Optimize prolong3: precompute coarse index/parity maps

(cherry picked from commit 2c0a3055d4)
This commit is contained in:
2026-02-28 23:53:30 +08:00
committed by ianchb
parent e7a02e8f72
commit 7109474a14

View File

@@ -1932,10 +1932,16 @@
integer,dimension(3) :: lbc,ubc,lbf,ubf,lbp,ubp,lbpc,ubpc integer,dimension(3) :: lbc,ubc,lbf,ubf,lbp,ubp,lbpc,ubpc
! when if=1 -> ic=0, this is different to vertex center grid ! when if=1 -> ic=0, this is different to vertex center grid
real*8, dimension(-2:extc(1),-2:extc(2),-2:extc(3)) :: funcc real*8, dimension(-2:extc(1),-2:extc(2),-2:extc(3)) :: funcc
integer,dimension(3) :: cxI integer,dimension(3) :: cxI
integer :: i,j,k,ii,jj,kk integer :: i,j,k,ii,jj,kk
real*8, dimension(6,6) :: tmp2 real*8, dimension(6,6) :: tmp2
real*8, dimension(6) :: tmp1 real*8, dimension(6) :: tmp1
integer, dimension(extf(1)) :: cix
integer, dimension(extf(2)) :: ciy
integer, dimension(extf(3)) :: ciz
logical, dimension(extf(1)) :: evenx
logical, dimension(extf(2)) :: eveny
logical, dimension(extf(3)) :: evenz
real*8, parameter :: C1=7.7d1/8.192d3,C2=-6.93d2/8.192d3,C3=3.465d3/4.096d3 real*8, parameter :: C1=7.7d1/8.192d3,C2=-6.93d2/8.192d3,C3=3.465d3/4.096d3
real*8, parameter :: C6=6.3d1/8.192d3,C5=-4.95d2/8.192d3,C4=1.155d3/4.096d3 real*8, parameter :: C6=6.3d1/8.192d3,C5=-4.95d2/8.192d3,C4=1.155d3/4.096d3
@@ -2002,10 +2008,10 @@
kmini=lbpc(3)-lbc(3) + 1 kmini=lbpc(3)-lbc(3) + 1
kmaxi=ubpc(3)-lbc(3) + 1 kmaxi=ubpc(3)-lbc(3) + 1
if(imino.lt.1.or.jmino.lt.1.or.kmino.lt.1.or.& if(imino.lt.1.or.jmino.lt.1.or.kmino.lt.1.or.&
imini.lt.1.or.jmini.lt.1.or.kmini.lt.1.or.& imini.lt.1.or.jmini.lt.1.or.kmini.lt.1.or.&
imaxo.gt.extf(1).or.jmaxo.gt.extf(2).or.kmaxo.gt.extf(3).or.& imaxo.gt.extf(1).or.jmaxo.gt.extf(2).or.kmaxo.gt.extf(3).or.&
imaxi.gt.extc(1)-2.or.jmaxi.gt.extc(2)-2.or.kmaxi.gt.extc(3)-2)then imaxi.gt.extc(1)-2.or.jmaxi.gt.extc(2)-2.or.kmaxi.gt.extc(3)-2)then
write(*,*)"error in prolongation for" write(*,*)"error in prolongation for"
write(*,*)"from" write(*,*)"from"
write(*,*)llbc,uubc write(*,*)llbc,uubc
@@ -2016,33 +2022,40 @@
write(*,*)"want" write(*,*)"want"
write(*,*)llbp,uubp write(*,*)llbp,uubp
write(*,*)lbp,ubp,lbpc,ubpc write(*,*)lbp,ubp,lbpc,ubpc
return return
endif endif
call symmetry_bd(3,extc,func,funcc,SoA) do i = imino,imaxo
ii = i + lbf(1) - 1
!~~~~~~> prolongation start... cix(i) = ii/2 - lbc(1) + 1
do k = kmino,kmaxo evenx(i) = (ii/2*2 == ii)
do j = jmino,jmaxo enddo
do i = imino,imaxo do j = jmino,jmaxo
cxI(1) = i jj = j + lbf(2) - 1
cxI(2) = j ciy(j) = jj/2 - lbc(2) + 1
cxI(3) = k eveny(j) = (jj/2*2 == jj)
! change to coarse level reference enddo
!|---*--- ---*--- ---*--- ---*--- ---*--- ---*--- ---*--- ---*---| do k = kmino,kmaxo
!|=======x===============x===============x===============x=======| kk = k + lbf(3) - 1
cxI = (cxI+lbf-1)/2 ciz(k) = kk/2 - lbc(3) + 1
! change to array index evenz(k) = (kk/2*2 == kk)
cxI = cxI - lbc + 1 enddo
if(any(cxI+3 > extc)) write(*,*)"error in prolong" call symmetry_bd(3,extc,func,funcc,SoA)
ii=i+lbf(1)-1
jj=j+lbf(2)-1 !~~~~~~> prolongation start...
kk=k+lbf(3)-1 do k = kmino,kmaxo
#if 0 do j = jmino,jmaxo
if(ii/2*2==ii)then do i = imino,imaxo
if(jj/2*2==jj)then cxI(1) = cix(i)
if(kk/2*2==kk)then cxI(2) = ciy(j)
cxI(3) = ciz(k)
if(any(cxI+3 > extc)) write(*,*)"error in prolong"
#if 0
if(ii/2*2==ii)then
if(jj/2*2==jj)then
if(kk/2*2==kk)then
tmp2= C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-2)+& tmp2= C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-2)+&
C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+& C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+&
C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+& C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+&
@@ -2125,11 +2138,11 @@
endif endif
endif endif
endif endif
#else #else
if(kk/2*2==kk)then if(evenz(k))then
tmp2= C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-2)+& tmp2= C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-2)+&
C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+& C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+&
C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+& C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+&
C4*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+1)+& C4*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+1)+&
C5*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+2)+& C5*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+2)+&
C6*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+3) C6*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+3)
@@ -2138,21 +2151,21 @@
C5*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+& C5*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)-1)+&
C4*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+& C4*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3) )+&
C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+1)+& C3*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+1)+&
C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+2)+& C2*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+2)+&
C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+3) C1*funcc(cxI(1)-2:cxI(1)+3,cxI(2)-2:cxI(2)+3,cxI(3)+3)
endif endif
if(jj/2*2==jj)then if(eveny(j))then
tmp1= C1*tmp2(:,1)+C2*tmp2(:,2)+C3*tmp2(:,3)+C4*tmp2(:,4)+C5*tmp2(:,5)+C6*tmp2(:,6) tmp1= C1*tmp2(:,1)+C2*tmp2(:,2)+C3*tmp2(:,3)+C4*tmp2(:,4)+C5*tmp2(:,5)+C6*tmp2(:,6)
else else
tmp1= C6*tmp2(:,1)+C5*tmp2(:,2)+C4*tmp2(:,3)+C3*tmp2(:,4)+C2*tmp2(:,5)+C1*tmp2(:,6) tmp1= C6*tmp2(:,1)+C5*tmp2(:,2)+C4*tmp2(:,3)+C3*tmp2(:,4)+C2*tmp2(:,5)+C1*tmp2(:,6)
endif endif
if(ii/2*2==ii)then if(evenx(i))then
funf(i,j,k)= C1*tmp1(1)+C2*tmp1(2)+C3*tmp1(3)+C4*tmp1(4)+C5*tmp1(5)+C6*tmp1(6) funf(i,j,k)= C1*tmp1(1)+C2*tmp1(2)+C3*tmp1(3)+C4*tmp1(4)+C5*tmp1(5)+C6*tmp1(6)
else else
funf(i,j,k)= C6*tmp1(1)+C5*tmp1(2)+C4*tmp1(3)+C3*tmp1(4)+C2*tmp1(5)+C1*tmp1(6) funf(i,j,k)= C6*tmp1(1)+C5*tmp1(2)+C4*tmp1(3)+C3*tmp1(4)+C2*tmp1(5)+C1*tmp1(6)
endif endif
#endif #endif
enddo enddo
enddo enddo