48 integer,
private :: itag
111 double precision,
dimension(:^D&,:),
allocatable :: w
121 integer :: nghostcellsCo, interpolation_order
122 integer :: nx^D, nxCo^D, ixG^L, i^D, idir
125 ixm^l=ixg^l^lsubnghostcells;
129 ixcogsmax^d=ixcogmax^d;
133 nx^d=ixmmax^d-ixmmin^d+1;
137 interpolation_order=1
139 interpolation_order=2
143 if (nghostcellsco+interpolation_order-1>
nghostcells)
then
144 call mpistop(
"interpolation order for prolongation in getbc too high")
150 ixs_srl_min^d(-1)=ixmmin^d
151 ixs_srl_min^d( 0)=ixmmin^d
154 ixs_srl_max^d( 0)=ixmmax^d
155 ixs_srl_max^d( 1)=ixmmax^d
158 ixr_srl_min^d( 0)=ixmmin^d
159 ixr_srl_min^d( 1)=ixmmax^d+1
161 ixr_srl_max^d( 0)=ixmmax^d
162 ixr_srl_max^d( 1)=ixgmax^d
164 ixs_r_min^d(-1)=ixcommin^d
165 ixs_r_min^d( 0)=ixcommin^d
168 ixs_r_max^d( 0)=ixcommax^d
169 ixs_r_max^d( 1)=ixcommax^d
172 ixr_r_min^d(1)=ixmmin^d
173 ixr_r_min^d(2)=ixmmin^d+nxco^d
174 ixr_r_min^d(3)=ixmmax^d+1
176 ixr_r_max^d(1)=ixmmin^d-1+nxco^d
177 ixr_r_max^d(2)=ixmmax^d
178 ixr_r_max^d(3)=ixgmax^d
180 ixs_p_min^d(0)=ixmmin^d-(interpolation_order-1)
181 ixs_p_min^d(1)=ixmmin^d-(interpolation_order-1)
182 ixs_p_min^d(2)=ixmmin^d+nxco^d-nghostcellsco-(interpolation_order-1)
183 ixs_p_min^d(3)=ixmmax^d+1-nghostcellsco-(interpolation_order-1)
184 ixs_p_max^d(0)=ixmmin^d-1+nghostcellsco+(interpolation_order-1)
185 ixs_p_max^d(1)=ixmmin^d-1+nxco^d+nghostcellsco+(interpolation_order-1)
186 ixs_p_max^d(2)=ixmmax^d+(interpolation_order-1)
187 ixs_p_max^d(3)=ixmmax^d+(interpolation_order-1)
191 ixs_p_min^d(0)=ixmmin^d
192 ixs_p_max^d(3)=ixmmax^d
193 ixs_p_max^d(1)=ixmmin^d-1+nxco^d+(interpolation_order-1)
194 ixs_p_min^d(2)=ixmmin^d+nxco^d-(interpolation_order-1)
197 ixr_p_min^d(0)=ixcommin^d-nghostcellsco-(interpolation_order-1)
198 ixr_p_min^d(1)=ixcommin^d-(interpolation_order-1)
199 ixr_p_min^d(2)=ixcommin^d-nghostcellsco-(interpolation_order-1)
200 ixr_p_min^d(3)=ixcommax^d+1-(interpolation_order-1)
202 ixr_p_max^d(1)=ixcommax^d+nghostcellsco+(interpolation_order-1)
203 ixr_p_max^d(2)=ixcommax^d+(interpolation_order-1)
204 ixr_p_max^d(3)=ixcommax^d+nghostcellsco+(interpolation_order-1)
208 ixr_p_min^d(3)=ixcommax^d+1
209 ixr_p_max^d(1)=ixcommax^d+(interpolation_order-1)
210 ixr_p_min^d(2)=ixcommin^d-(interpolation_order-1)
216 allocate(pole_buf%ws(ixgs^t,nws))
219 { ixs_srl_stg_min^d(idir,-1)=ixmmin^d-
kr(idir,^d)
221 ixs_srl_stg_min^d(idir,0) =ixmmin^d-
kr(idir,^d)
222 ixs_srl_stg_max^d(idir,0) =ixmmax^d
223 ixs_srl_stg_min^d(idir,1) =ixmmax^d-
nghostcells+1-
kr(idir,^d)
224 ixs_srl_stg_max^d(idir,1) =ixmmax^d
226 ixr_srl_stg_min^d(idir,-1)=1-
kr(idir,^d)
228 ixr_srl_stg_min^d(idir,0) =ixmmin^d-
kr(idir,^d)
229 ixr_srl_stg_max^d(idir,0) =ixmmax^d
230 ixr_srl_stg_min^d(idir,1) =ixmmax^d+1-
kr(idir,^d)
231 ixr_srl_stg_max^d(idir,1) =ixgmax^d
233 ixs_r_stg_min^d(idir,-1)=ixcommin^d-
kr(idir,^d)
235 ixs_r_stg_min^d(idir,0) =ixcommin^d-
kr(idir,^d)
236 ixs_r_stg_max^d(idir,0) =ixcommax^d
237 ixs_r_stg_min^d(idir,1) =ixcommax^d+1-
nghostcells-
kr(idir,^d)
238 ixs_r_stg_max^d(idir,1) =ixcommax^d
240 ixr_r_stg_min^d(idir,0)=1-
kr(idir,^d)
242 ixr_r_stg_min^d(idir,1)=ixmmin^d-
kr(idir,^d)
243 ixr_r_stg_max^d(idir,1)=ixmmin^d-1+nxco^d
244 ixr_r_stg_min^d(idir,2)=ixmmin^d+nxco^d-
kr(idir,^d)
245 ixr_r_stg_max^d(idir,2)=ixmmax^d
246 ixr_r_stg_min^d(idir,3)=ixmmax^d+1-
kr(idir,^d)
247 ixr_r_stg_max^d(idir,3)=ixgmax^d
252 ixs_p_stg_min^d(idir,0)=ixmmin^d-1
253 ixs_p_stg_max^d(idir,0)=ixmmin^d-1+nghostcellsco
254 ixs_p_stg_min^d(idir,1)=ixmmin^d-1
255 ixs_p_stg_max^d(idir,1)=ixmmin^d-1+nxco^d+nghostcellsco
256 ixs_p_stg_min^d(idir,2)=ixmmax^d-nxco^d-nghostcellsco
257 ixs_p_stg_max^d(idir,2)=ixmmax^d
258 ixs_p_stg_min^d(idir,3)=ixmmax^d-nghostcellsco
259 ixs_p_stg_max^d(idir,3)=ixmmax^d
261 ixr_p_stg_min^d(idir,0)=ixcommin^d-1-nghostcellsco
262 ixr_p_stg_max^d(idir,0)=ixcommin^d-1
263 ixr_p_stg_min^d(idir,1)=ixcommin^d-1
264 ixr_p_stg_max^d(idir,1)=ixcommax^d+nghostcellsco
265 ixr_p_stg_min^d(idir,2)=ixcommin^d-1-nghostcellsco
266 ixr_p_stg_max^d(idir,2)=ixcommax^d
267 ixr_p_stg_min^d(idir,3)=ixcommax^d+1-1
268 ixr_p_stg_max^d(idir,3)=ixcommax^d+nghostcellsco
273 ixs_p_stg_min^d(idir,0)=ixmmin^d
274 ixs_p_stg_max^d(idir,0)=ixmmin^d-1+nghostcellsco+(interpolation_order-1)
275 ixs_p_stg_min^d(idir,1)=ixmmin^d
276 ixs_p_stg_max^d(idir,1)=ixmmin^d-1+nxco^d+nghostcellsco+(interpolation_order-1)
277 ixs_p_stg_min^d(idir,2)=ixmmax^d+1-nxco^d-nghostcellsco-(interpolation_order-1)
278 ixs_p_stg_max^d(idir,2)=ixmmax^d
279 ixs_p_stg_min^d(idir,3)=ixmmax^d+1-nghostcellsco-(interpolation_order-1)
280 ixs_p_stg_max^d(idir,3)=ixmmax^d
282 ixr_p_stg_min^d(idir,0)=ixcommin^d-nghostcellsco-(interpolation_order-1)
283 ixr_p_stg_max^d(idir,0)=ixcommin^d-1
284 ixr_p_stg_min^d(idir,1)=ixcommin^d
285 ixr_p_stg_max^d(idir,1)=ixcommax^d+nghostcellsco+(interpolation_order-1)
286 ixr_p_stg_min^d(idir,2)=ixcommin^d-nghostcellsco-(interpolation_order-1)
287 ixr_p_stg_max^d(idir,2)=ixcommax^d
288 ixr_p_stg_min^d(idir,3)=ixcommax^d+1
289 ixr_p_stg_max^d(idir,3)=ixcommax^d+nghostcellsco+(interpolation_order-1)
298 sizes_srl_send_stg(idir,i^d)={(ixs_srl_stg_max^d(idir,i^d)-ixs_srl_stg_min^d(idir,i^d)+1)|*}
299 sizes_srl_recv_stg(idir,i^d)={(ixr_srl_stg_max^d(idir,i^d)-ixr_srl_stg_min^d(idir,i^d)+1)|*}
300 sizes_r_send_stg(idir,i^d)={(ixs_r_stg_max^d(idir,i^d)-ixs_r_stg_min^d(idir,i^d)+1)|*}
310 sizes_r_recv_stg(idir,i^d)={(ixr_r_stg_max^d(idir,i^d)-ixr_r_stg_min^d(idir,i^d)+1)|*}
311 sizes_p_send_stg(idir,i^d)={(ixs_p_stg_max^d(idir,i^d)-ixs_p_stg_min^d(idir,i^d)+1)|*}
312 sizes_p_recv_stg(idir,i^d)={(ixr_p_stg_max^d(idir,i^d)-ixr_p_stg_min^d(idir,i^d)+1)|*}
325 integer,
intent(in) :: nwstart, nwbc
326 integer :: i^D, ic^D, inc^D
329 if(i^d==0|.and.) cycle
333 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
334 inc^db=2*i^db+ic^db\}
346 integer,
intent(inout) :: comm_type
347 integer,
intent(in) :: ix^L, ixG^L, nwstart, nwbc
349 integer,
dimension(ndim+1) :: fullsize, subsize, start
351 ^
d&fullsize(^
d)=ixgmax^
d;
353 ^
d&subsize(^
d)=ixmax^
d-ixmin^
d+1;
355 ^
d&start(^
d)=ixmin^
d-1;
356 start(
ndim+1)=nwstart-1
358 call mpi_type_create_subarray(
ndim+1,fullsize,subsize,start,mpi_order_fortran, &
359 mpi_double_precision,comm_type,
ierrmpi)
360 call mpi_type_commit(comm_type,
ierrmpi)
365 subroutine getbc(time,qdt,psb,nwstart,nwbc,req_diag)
372 double precision,
intent(in) :: time, qdt
373 type(state),
target :: psb(max_blocks)
374 integer,
intent(in) :: nwstart
375 integer,
intent(in) :: nwbc
376 logical,
intent(in),
optional :: req_diag
378 double precision :: time_bcin
379 integer :: ipole, nwhead, nwtail
380 integer :: iigrid, igrid, ineighbor, ipe_neighbor, isizes
381 integer :: ixR^L, ixS^L
382 integer :: i^D, n_i^D, ic^D, inc^D, n_inc^D, idir
383 integer :: isend_buf(npwbuf), ipwbuf, nghostcellsco
385 integer :: ibuf_start, ibuf_next
387 integer,
dimension(1) :: shapes
388 logical :: req_diagonal
391 time_bcin=mpi_wtime()
394 nwtail=nwstart+nwbc-1
396 req_diagonal = .true.
397 if (
present(req_diag)) req_diagonal = req_diag
406 do iigrid=1,igridstail; igrid=igrids(iigrid);
407 if(any(neighbor_type(:^d&,igrid)==neighbor_coarse))
then
436 do iigrid=1,igridstail; igrid=igrids(iigrid);
439 select case (neighbor_type(i^d,igrid))
440 case (neighbor_sibling)
449 do iigrid=1,igridstail; igrid=igrids(iigrid);
452 select case (neighbor_type(i^d,igrid))
453 case (neighbor_sibling)
455 case (neighbor_coarse)
463 do iigrid=1,igridstail; igrid=igrids(iigrid);
466 select case (neighbor_type(i^d,igrid))
467 case(neighbor_sibling)
469 case(neighbor_coarse)
479 if(stagger_grid)
then
487 do iigrid=1,igridstail; igrid=igrids(iigrid);
490 select case (neighbor_type(i^d,igrid))
491 case (neighbor_sibling)
501 if (isend_buf(ipwbuf)/=0)
deallocate(pwbuf(ipwbuf)%w)
510 do iigrid=1,igridstail; igrid=igrids(iigrid);
517 do iigrid=1,igridstail; igrid=igrids(iigrid);
526 do iigrid=1,igridstail; igrid=igrids(iigrid);
529 if (neighbor_type(i^d,igrid)==neighbor_fine)
call bc_fill_prolong(igrid,i^d)
537 if(stagger_grid)
then
543 do iigrid=1,igridstail; igrid=igrids(iigrid);
552 do iigrid=1,igridstail; igrid=igrids(iigrid);
558 if (isend_buf(ipwbuf)/=0)
deallocate(pwbuf(ipwbuf)%w)
564 do iigrid=1,igridstail; igrid=igrids(iigrid);
565 if(.not.phyboundblock(igrid)) cycle
572 if(
bcphys.and.
associated(phys_boundary_adjust))
then
574 do iigrid=1,igridstail; igrid=igrids(iigrid);
575 if(.not.phyboundblock(igrid)) cycle
576 call phys_boundary_adjust(igrid,psb)
581 time_bc=time_bc+(mpi_wtime()-time_bcin)
586 integer,
intent(in) :: dir(^nd)
588 if (all(dir == 0))
then
590 else if (.not. req_diagonal .and. count(dir /= 0) > 1)
then
600 integer,
intent(in) :: igrid
602 integer :: idims,iside,i^D,k^L,ixB^L
605 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
612 kmin2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0,-1,igrid)==1)
613 kmax2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0, 1,igrid)==1)}
615 kmin2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0,-1,0,igrid)==1)
616 kmax2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0, 1,0,igrid)==1)
617 kmin3=merge(1, 0, idims .lt. 3 .and. neighbor_type(0,0,-1,igrid)==1)
618 kmax3=merge(1, 0, idims .lt. 3 .and. neighbor_type(0,0, 1,igrid)==1)}
619 ixbmin^d=ixglo^d+kmin^d*nghostcells;
620 ixbmax^d=ixghi^d-kmax^d*nghostcells;
622 i^d=kr(^d,idims)*(2*iside-3);
623 if (aperiodb(idims))
then
624 if (neighbor_type(i^d,igrid) /= neighbor_boundary .and. &
625 .not. psb(igrid)%is_physical_boundary(2*idims-2+iside)) cycle
627 if (neighbor_type(i^d,igrid) /= neighbor_boundary) cycle
629 call bc_phys(iside,idims,time,qdt,psb(igrid),ixg^ll,ixb^l)
638 ipe_neighbor=neighbor(2,i^d,igrid)
639 if (ipe_neighbor/=mype)
then
641 itag=(3**^nd+4**^nd)*(igrid-1)+{(i^d+1)*3**(^d-1)+}
644 if(stagger_grid)
then
657 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
658 inc^db=2*i^db+ic^db\}
659 ipe_neighbor=neighbor_child(2,inc^d,igrid)
660 if (ipe_neighbor/=mype)
then
662 itag=(3**^nd+4**^nd)*(igrid-1)+3**^nd+{inc^d*4**(^d-1)+}
663 call mpi_irecv(psb(igrid)%w,1,
type_recv_r(inc^d), &
665 if(stagger_grid)
then
668 mpi_double_precision,ipe_neighbor,itag, &
680 ipe_neighbor=neighbor(2,i^d,igrid)
682 if(ipe_neighbor/=mype)
then
683 ineighbor=neighbor(1,i^d,igrid)
684 ipole=neighbor_pole(i^d,igrid)
688 itag=(3**^nd+4**^nd)*(ineighbor-1)+{(n_i^d+1)*3**(^d-1)+}
691 if(stagger_grid)
then
698 reshape(psb(igrid)%ws(ixs^s,idir),shapes)
710 n_i^d=i^d^d%n_i^dd=-i^dd;\}
712 if (isend_buf(ipwbuf)/=0)
then
715 deallocate(pwbuf(ipwbuf)%w)
717 allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
718 call pole_buffer(pwbuf(ipwbuf)%w,ixs^
l,ixs^
l,psb(igrid)%w,ixg^ll,ixs^
l)
721 itag=(3**^nd+4**^nd)*(ineighbor-1)+{(n_i^d+1)*3**(^d-1)+}
722 isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
723 call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
725 ipwbuf=1+modulo(ipwbuf,
npwbuf)
726 if(stagger_grid)
then
733 reshape(psb(igrid)%ws(ixs^s,idir),shapes)
747 integer,
intent(in) :: igrid,i^D
748 integer :: ineighbor,ipe_neighbor,ipole,ixS^L,ixR^L,n_i^D,idir
750 ipe_neighbor=neighbor(2,i^d,igrid)
751 if(ipe_neighbor==mype)
then
752 ineighbor=neighbor(1,i^d,igrid)
753 ipole=neighbor_pole(i^d,igrid)
758 psb(ineighbor)%w(ixr^s,nwhead:nwtail)=&
759 psb(igrid)%w(ixs^s,nwhead:nwtail)
760 if(stagger_grid)
then
764 psb(ineighbor)%ws(ixr^s,idir)=psb(igrid)%ws(ixs^s,idir)
771 n_i^d=i^d^d%n_i^dd=-i^dd;\}
774 call pole_copy(psb(ineighbor)%w,ixg^ll,ixr^l,psb(igrid)%w,ixg^ll,ixs^l,ipole)
775 if(stagger_grid)
then
779 call pole_copy_stg(psb(ineighbor)%ws,ixgs^ll,ixr^l,psb(igrid)%ws,ixgs^ll,ixs^l,idir,ipole)
790 ipe_neighbor=neighbor(2,i^d,igrid)
791 if(ipe_neighbor/=mype)
then
792 ic^d=1+modulo(node(pig^d_,igrid)-1,2);
793 if({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.})
return
794 ineighbor=neighbor(1,i^d,igrid)
795 ipole=neighbor_pole(i^d,igrid)
799 itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
802 if(stagger_grid)
then
809 reshape(psc(igrid)%ws(ixs^s,idir),shapes)
814 mpi_double_precision,ipe_neighbor,itag, &
822 n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
824 if(isend_buf(ipwbuf)/=0)
then
827 deallocate(pwbuf(ipwbuf)%w)
829 allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
833 itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
834 isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
835 call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
837 ipwbuf=1+modulo(ipwbuf,
npwbuf)
838 if(stagger_grid)
then
845 reshape(psc(igrid)%ws(ixs^s,idir),shapes)
850 mpi_double_precision,ipe_neighbor,itag, &
861 integer,
intent(in) :: igrid,i^D
863 integer :: ic^D,n_inc^D,ixS^L,ixR^L,ipe_neighbor,ineighbor,ipole,idir
865 ipe_neighbor=neighbor(2,i^d,igrid)
866 if(ipe_neighbor==mype)
then
867 ic^d=1+modulo(node(pig^d_,igrid)-1,2);
868 if({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.})
return
869 ineighbor=neighbor(1,i^d,igrid)
870 ipole=neighbor_pole(i^d,igrid)
875 psb(ineighbor)%w(ixr^s,nwhead:nwtail)=&
876 psc(igrid)%w(ixs^s,nwhead:nwtail)
877 if(stagger_grid)
then
881 psb(ineighbor)%ws(ixr^s,idir)=psc(igrid)%ws(ixs^s,idir)
888 n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
891 call pole_copy(psb(ineighbor)%w,ixg^ll,ixr^l,psc(igrid)%w,
ixcog^l,ixs^l,ipole)
892 if(stagger_grid)
then
897 call pole_copy_stg(psb(ineighbor)%ws,ixgs^ll,ixr^l,psc(igrid)%ws,
ixcogs^l,ixs^l,idir,ipole)
907 integer :: ixS^L,ixR^L,n_i^D,ixSsync^L,ixRsync^L,idir
909 ipe_neighbor=neighbor(2,i^d,igrid)
910 if(ipe_neighbor/=mype)
then
911 ineighbor=neighbor(1,i^d,igrid)
912 ipole=neighbor_pole(i^d,igrid)
924 shape=shape(psb(igrid)%ws(ixs^s,idir)))
930 n_i^d=i^d^d%n_i^dd=-i^dd;\}
938 shape=shape(psb(igrid)%ws(ixs^s,idir)))
940 call pole_copy_stg(psb(igrid)%ws,ixgs^ll,ixr^l,pole_buf%ws,ixgs^ll,ixs^l,idir,ipole)
950 ipole=neighbor_pole(i^d,igrid)
953 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
954 inc^db=2*i^db+ic^db\}
955 ipe_neighbor=neighbor_child(2,inc^d,igrid)
956 if(ipe_neighbor/=mype)
then
957 ineighbor=neighbor_child(1,inc^d,igrid)
964 shape=shape(psb(igrid)%ws(ixr^s,idir)))
970 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
971 inc^db=2*i^db+ic^db\}
972 ipe_neighbor=neighbor_child(2,inc^d,igrid)
973 if(ipe_neighbor/=mype)
then
974 ineighbor=neighbor_child(1,inc^d,igrid)
977 n_i^d=i^d^d%n_i^dd=-i^dd;\}
987 shape=shape(psb(igrid)%ws(ixr^s,idir)))
988 call pole_copy_stg(psb(igrid)%ws,ixgs^ll,ixr^
l,pole_buf%ws,ixgs^ll,ixr^
l,idir,ipole)
1000 ic^d=1+modulo(node(pig^d_,igrid)-1,2);
1001 if ({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.})
return
1003 ipe_neighbor=neighbor(2,i^d,igrid)
1004 if (ipe_neighbor/=mype)
then
1007 itag=(3**^nd+4**^nd)*(igrid-1)+3**^nd+{inc^d*4**(^d-1)+}
1008 call mpi_irecv(psc(igrid)%w,1,
type_recv_p(inc^d), &
1010 if(stagger_grid)
then
1013 mpi_double_precision,ipe_neighbor,itag,&
1025 ipole=neighbor_pole(i^d,igrid)
1027 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1028 inc^db=2*i^db+ic^db\}
1029 ipe_neighbor=neighbor_child(2,inc^d,igrid)
1030 if(ipe_neighbor/=mype)
then
1031 ineighbor=neighbor_child(1,inc^d,igrid)
1036 itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
1037 call mpi_isend(psb(igrid)%w,1,
type_send_p(inc^d), &
1039 if(stagger_grid)
then
1046 reshape(psb(igrid)%ws(ixs^s,idir),shapes)
1047 ibuf_start=ibuf_next
1051 mpi_double_precision,ipe_neighbor,itag, &
1059 n_inc^d=inc^d^d%n_inc^dd=ic^dd-i^dd;\}
1061 if(isend_buf(ipwbuf)/=0)
then
1064 deallocate(pwbuf(ipwbuf)%w)
1066 allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
1067 call pole_buffer(pwbuf(ipwbuf)%w,ixs^
l,ixs^
l,psb(igrid)%w,ixg^ll,ixs^
l)
1070 itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
1071 isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
1072 call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
1074 ipwbuf=1+modulo(ipwbuf,
npwbuf)
1075 if(stagger_grid)
then
1082 reshape(psb(igrid)%ws(ixs^s,idir),shapes)
1083 ibuf_start=ibuf_next
1087 mpi_double_precision,ipe_neighbor,itag, &
1099 integer,
intent(in) :: igrid,i^D
1101 integer :: ipe_neighbor,ineighbor,ixS^L,ixR^L,ic^D,inc^D,ipole,idir
1103 ipole=neighbor_pole(i^d,igrid)
1106 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1107 inc^db=2*i^db+ic^db\}
1108 ipe_neighbor=neighbor_child(2,inc^d,igrid)
1109 if(ipe_neighbor==mype)
then
1111 ineighbor=neighbor_child(1,inc^d,igrid)
1115 psc(ineighbor)%w(ixr^s,nwhead:nwtail) &
1116 =psb(igrid)%w(ixs^s,nwhead:nwtail)
1117 if(stagger_grid)
then
1121 psc(ineighbor)%ws(ixr^s,idir)=psb(igrid)%ws(ixs^s,idir)
1127 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1128 inc^db=2*i^db+ic^db\}
1129 ipe_neighbor=neighbor_child(2,inc^d,igrid)
1130 if(ipe_neighbor==mype)
then
1132 ineighbor=neighbor_child(1,inc^d,igrid)
1135 n_inc^d=inc^d^d%n_inc^dd=ic^dd-i^dd;\}
1139 if(stagger_grid)
then
1152 integer,
intent(in) :: igrid
1154 integer :: i^D,idims,iside
1155 logical,
dimension(-1:1^D&) :: NeedProlong
1160 if (neighbor_type(i^d,igrid)==neighbor_coarse)
then
1162 needprolong(i^d)=.true.
1165 if(stagger_grid)
then
1176 if (needprolong(i^dd))
call bc_prolong_stg(igrid,i^dd,needprolong)
1212 ic^d=1+modulo(node(pig^d_,igrid)-1,2);
1213 if ({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.})
return
1215 ipe_neighbor=neighbor(2,i^d,igrid)
1216 if(ipe_neighbor/=mype)
then
1217 ineighbor=neighbor(1,i^d,igrid)
1218 ipole=neighbor_pole(i^d,igrid)
1227 shape=shape(psc(igrid)%ws(ixr^s,idir)))
1234 n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
1242 shape=shape(psc(igrid)%ws(ixr^s,idir)))
1255 double precision :: dxFi^D, dxCo^D, xFimin^D, xComin^D, invdxCo^D
1256 integer :: i^D,igrid
1257 integer :: ixFi^L,ixCo^L,ii^D, idims,iside,ixB^L
1260 dxfi^d=rnode(rpdx^d_,igrid);
1262 invdxco^d=1.d0/dxco^d;
1268 xfimin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxfi^d;
1269 xcomin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxco^d;
1271 if(phyboundblock(igrid).and.
bcphys)
then
1274 ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1275 ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1279 if(neighbor_type(-1,0,0,igrid)==neighbor_boundary .or. &
1280 neighbor_type(1,0,0,igrid)==neighbor_boundary)
then
1281 if(neighbor_type(0,-1,0,igrid)==neighbor_boundary) ixcomin2=ixcommin2
1282 if(neighbor_type(0,0,-1,igrid)==neighbor_boundary) ixcomin3=ixcommin3
1283 if(neighbor_type(0,1,0,igrid)==neighbor_boundary) ixcomax2=ixcommax2
1284 if(neighbor_type(0,0,1,igrid)==neighbor_boundary) ixcomax3=ixcommax3
1286 else if(idims == 2)
then
1287 if(neighbor_type(0,-1,0,igrid)==neighbor_boundary .or. &
1288 neighbor_type(0,1,0,igrid)==neighbor_boundary)
then
1289 if(neighbor_type(0,0,-1,igrid)==neighbor_boundary) ixcomin3=ixcommin3
1290 if(neighbor_type(0,0,1,igrid)==neighbor_boundary) ixcomax3=ixcommax3
1295 ii^d=kr(^d,idims)*(2*iside-3);
1296 if(neighbor_type(ii^d,igrid)/=neighbor_boundary) cycle
1297 if(( {(iside==1.and.idims==^d.and.ixcomin^d<ixcogmin^d+nghostcells)|.or.} ) &
1298 .or.( {(iside==2.and.idims==^d.and.ixcomax^d>ixcogmax^d-nghostcells)|.or. }))
then
1299 {ixbmin^d=merge(ixcogmin^d,ixcomin^d,idims==^d);}
1300 {ixbmax^d=merge(ixcogmax^d,ixcomax^d,idims==^d);}
1301 call bc_phys(iside,idims,time,0.d0,psc(igrid),
ixcog^l,ixb^l)
1307 if(prolongprimitive)
then
1314 ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1315 ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1325 if(prolongprimitive)
then
1334 double precision :: dxFi^D,dxCo^D,xFimin^D,xComin^D,invdxCo^D
1335 integer :: igrid,i^D
1336 integer :: ixFi^L,ixCo^L
1337 logical,
dimension(-1:1^D&) :: NeedProlong
1338 logical :: fine_^Lin
1342 if(i^d>-1) fine_min^din=(.not.needprolong(i^dd-kr(^d,^dd)).and.neighbor_type(i^dd-kr(^d,^dd),igrid)/=1)
1343 if(i^d<1) fine_max^din=(.not.needprolong(i^dd+kr(^d,^dd)).and.neighbor_type(i^dd+kr(^d,^dd),igrid)/=1)
1348 dxfi^d=rnode(rpdx^d_,igrid);
1350 invdxco^d=1.d0/dxco^d;
1352 xfimin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxfi^d;
1353 xcomin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxco^d;
1358 ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1359 ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1361 if(prolongprimitive)
call phys_to_primitive(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1363 call prolong_2nd_stg(psc(igrid),psb(igrid),ixco^l,ixfi^l,dxco^d,xcomin^d,dxfi^d,xfimin^d,.true.,fine_^lin)
1365 if(prolongprimitive)
call phys_to_conserved(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1368 needprolong(i^d)=.false.
1373 dxCo^D,invdxCo^D,xComin^D)
1375 integer,
intent(in) :: igrid, ixFi^L
1376 double precision,
intent(in) :: dxFi^D, xFimin^D,dxCo^D, invdxCo^D, xComin^D
1378 double precision :: xCo^D, xFi^D, eta^D
1379 double precision :: slopeL, slopeR, slopeC, signC, signR
1380 double precision :: slope(1:nw,ndim)
1382 double precision :: signedfactorhalf^D
1383 integer :: ixCo^D, jxCo^D, hxCo^D, ixFi^D, ix^D, iw, idims, nwmin,nwmax
1388 if(prolongprimitive)
then
1396 {
do ixfi^db = ixfi^lim^db
1399 xfi^db=xfimin^db+(dble(ixfi^db)-half)*dxfi^db
1404 ixco^db=int((xfi^db-xcomin^db)*invdxco^db)+1
1408 xco^db=xcomin^db+(dble(ixco^db)-half)*dxco^db \}
1414 if(slab_uniform)
then
1423 eta^d=(xfi^d-xco^d)*invdxco^d;
1457 ix^d=2*int((ixfi^d+ixmlo^d)/2)-ixmlo^d;
1458 {
if(xfi^d>xco^d)
then
1459 signedfactorhalf^d=0.5d0
1461 signedfactorhalf^d=-0.5d0
1463 eta^d=signedfactorhalf^d*(one-psb(igrid)%dvolume(ixfi^dd) &
1464 /sum(psb(igrid)%dvolume(ix^d:ix^d+1^d%ixFi^dd))) \}
1471 hxco^d=ixco^d-kr(^d,idims)\
1472 jxco^d=ixco^d+kr(^d,idims)\
1475 slopel=psc(igrid)%w(ixco^d,iw)-psc(igrid)%w(hxco^d,iw)
1476 sloper=psc(igrid)%w(jxco^d,iw)-psc(igrid)%w(ixco^d,iw)
1477 slopec=half*(sloper+slopel)
1480 signr=sign(one,sloper)
1481 signc=sign(one,slopec)
1499 slope(iw,idims)=signc*max(zero,min(dabs(slopec), &
1500 signc*slopel,signc*sloper))
1506 psb(igrid)%w(ixfi^d,nwmin:nwmax)=psc(igrid)%w(ixco^d,nwmin:nwmax)+&
1507 {(slope(nwmin:nwmax,^d)*eta^d)+}
1511 if(prolongprimitive)
then
1513 call phys_to_conserved(ixg^ll,ixfi^
l,psb(igrid)%w,psb(igrid)%x)
1519 dxCo^D,invdxCo^D,xComin^D)
1521 integer,
intent(in) :: igrid, ixFi^L
1522 double precision,
intent(in) :: dxFi^D, xFimin^D,dxCo^D, invdxCo^D, xComin^D
1524 double precision :: xFi^D
1525 integer :: ixCo^D, ixFi^D, nwmin,nwmax
1527 if(prolongprimitive)
then
1535 {
do ixfi^db = ixfi^lim^db
1537 xfi^db=xfimin^db+(dble(ixfi^db)-half)*dxfi^db
1541 ixco^db=int((xfi^db-xcomin^db)*invdxco^db)+1\}
1544 psb(igrid)%w(ixfi^d,nwmin:nwmax)=psc(igrid)%w(ixco^d,nwmin:nwmax)
1548 if(prolongprimitive)
call phys_to_conserved(ixg^ll,ixfi^
l,psb(igrid)%w,psb(igrid)%x)
1552 subroutine pole_copy(wrecv,ixIR^L,ixR^L,wsend,ixIS^L,ixS^L,ipole)
1554 integer,
intent(in) :: ixIR^L,ixR^L,ixIS^L,ixS^L,ipole
1555 double precision :: wrecv(ixIR^S,1:nw), wsend(ixIS^S,1:nw)
1557 integer :: iw, iside, iB
1561 iside=int((i^d+3)/2)
1564 select case (typeboundary(iw,ib))
1566 wrecv(ixr^s,iw) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1568 wrecv(ixr^s,iw) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1570 call mpistop(
"Pole boundary condition should be symm or asymm")
1579 integer,
intent(in) :: ixIR^L,ixR^L,ixIS^L,ixS^L,idirs,ipole
1581 double precision :: wrecv(ixIR^S,1:nws), wsend(ixIS^S,1:nws)
1582 integer :: iB, iside
1586 iside=int((i^d+3)/2)
1588 select case (typeboundary(iw_mag(idirs),ib))
1590 wrecv(ixr^s,idirs) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,idirs)
1592 wrecv(ixr^s,idirs) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,idirs)
1594 call mpistop(
"Pole boundary condition should be symm or asymm")
1603 integer,
intent(in) :: ixIR^L,ixR^L,ixIS^L,ixS^L
1604 double precision :: wrecv(ixIR^S,nwhead:nwtail), wsend(ixIS^S,1:nw)
1606 integer :: iw, iside, iB
1610 iside=int((i^d+3)/2)
1613 select case (typeboundary(iw,ib))
1615 wrecv(ixr^s,iw) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1617 wrecv(ixr^s,iw) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1619 call mpistop(
"Pole boundary condition should be symm or asymm")
1626 end subroutine getbc
subroutine bc_recv_restrict
Receive from fine neighbor.
subroutine gc_prolong(igrid)
subroutine interpolation_linear(igrid, ixFiL, dxFiD, xFiminD, dxCoD, invdxCoD, xCominD)
subroutine bc_fill_srl(igrid, iD)
subroutine bc_fill_srl_stg
fill siblings ghost cells with received data
subroutine pole_copy(wrecv, ixIRL, ixRL, wsend, ixISL, ixSL, ipole)
subroutine bc_recv_prolong
Receive from coarse neighbor.
subroutine bc_fill_prolong(igrid, iD)
Send to finer neighbor.
subroutine bc_prolong_stg(igrid, iD, NeedProlong)
subroutine bc_recv_srl
Receive from sibling at same refinement level.
subroutine interpolation_copy(igrid, ixFiL, dxFiD, xFiminD, dxCoD, invdxCoD, xCominD)
subroutine pole_copy_stg(wrecv, ixIRL, ixRL, wsend, ixISL, ixSL, idirs, ipole)
subroutine bc_prolong(igrid, iD)
do prolongation for fine blocks after receipt data from coarse neighbors
subroutine bc_send_restrict
Send to coarser neighbor.
subroutine pole_buffer(wrecv, ixIRL, ixRL, wsend, ixISL, ixSL)
subroutine bc_fill_restrict(igrid, iD)
fill coarser neighbor's ghost cells
subroutine bc_fill_restrict_stg
fill restricted ghost cells after receipt
subroutine bc_send_srl
Send to sibling at same refinement level.
subroutine fill_boundary_after_gc(igrid)
Physical boundary conditions.
logical function skip_direction(dir)
subroutine bc_send_prolong
Send to finer neighbor.
subroutine bc_fill_prolong_stg
fill coarser representative with data from coarser neighbors
subroutine, public prolong_2nd_stg(sCo, sFi, ixCoLin, ixFiLin, dxCoD, xCominD, dxFiD, xFiminD, ghost, fine_Lin)
This subroutine performs a 2nd order prolongation for a staggered field F, preserving the divergence ...
subroutine, public getintbc(time, ixGL)
fill inner boundary values
subroutine, public bc_phys(iside, idims, time, qdt, s, ixGL, ixBL)
fill ghost cells at a physical boundary
subroutine, public coarsen_grid(sFi, ixFiGL, ixFiL, sCo, ixCoGL, ixCoL)
coarsen one grid to its coarser representative
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
update ghost cells of all blocks including physical boundaries
integer, dimension(0:3^d &), target type_recv_p_p1
integer, dimension( :^d &), pointer type_recv_r
integer, dimension(-1:1^d &), target type_send_srl_p1
integer, dimension(-1:1^d &) sizes_r_send_total
integer, dimension(^nd, 0:3) l
integer, dimension(-1:1^d &), target type_send_r_p2
integer, dimension(:), allocatable sendrequest_c_p
integer, dimension(^nd,-1:1) ixs_r_stg_
integer, dimension(^nd,-1:1) ixs_srl_stg_
subroutine get_bc_comm_type(comm_type, ixL, ixGL, nwstart, nwbc)
integer, dimension(^nd, 0:3^d &) sizes_p_send_stg
integer, dimension(-1:1^d &) sizes_srl_send_total
integer, dimension(^nd, 0:3) ixr_r_stg_
integer, dimension(:), allocatable sendrequest_r
integer, dimension(:), allocatable recvrequest_r
double precision, dimension(:), allocatable sendbuffer_p
integer, dimension(:,:), allocatable sendstatus_c_p
integer, dimension(^nd,-1:1^d &) sizes_r_send_stg
integer, dimension(0:3) ixs_p_
integer, dimension(0:3) ixr_r_
integer, dimension(0:3^d &), target type_send_p_f
integer, dimension( :^d &), pointer type_send_p
integer, dimension(:,:), allocatable sendstatus_c_sr
integer, dimension(0:3^d &), target type_recv_p_p2
integer, dimension(0:3) ixr_p_
double precision, dimension(:), allocatable recvbuffer_r
integer, dimension( :^d &), pointer type_send_srl
integer, dimension(-1:1^d &), target type_send_r_p1
integer, dimension(:), allocatable recvrequest_c_p
integer, dimension(-1:1^d &), target type_send_srl_p2
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(0:3^d &), target type_send_p_p2
integer, dimension(:), allocatable recvrequest_c_sr
integer, dimension(0:3^d &), target type_recv_r_f
integer, dimension(-1:1^d &), target type_recv_srl_f
integer, dimension(:), allocatable sendrequest_c_sr
integer, dimension(-1:1^d &), target type_send_r_f
double precision, dimension(:), allocatable recvbuffer_srl
integer, dimension(:,:), allocatable recvstatus_p
integer, dimension(^nd, 0:3) ixs_p_stg_
integer, dimension(-1:1^d &), target type_send_srl_f
integer, dimension(^nd, 0:3^d &) sizes_p_recv_stg
integer, dimension(0:3^d &), target type_send_p_p1
integer, dimension(^nd,-1:1^d &) sizes_srl_recv_stg
double precision, dimension(:), allocatable recvbuffer_p
integer, dimension(^nd,-1:1) ixr_srl_stg_
integer, dimension(-1:1) ixs_r_
integer, dimension( :^d &), pointer type_send_r
integer, dimension(-1:1) ixr_srl_
integer, dimension(0:3^d &), target type_recv_r_p1
double precision, dimension(:), allocatable sendbuffer_r
integer, dimension(0:3^d &), target type_recv_r_p2
integer, dimension(0:3^d &) sizes_p_recv_total
integer, dimension(:,:), allocatable recvstatus_c_sr
integer, dimension(0:3^d &) sizes_r_recv_total
integer, parameter npwbuf
integer, dimension(:), allocatable sendrequest_p
integer, dimension(^nd, 0:3^d &) sizes_r_recv_stg
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
integer, dimension(-1:1^d &), target type_recv_srl_p2
double precision, dimension(:), allocatable sendbuffer_srl
integer, dimension(-1:1) ixs_srl_
integer, dimension(:,:), allocatable sendstatus_srl
integer, dimension( :^d &), pointer type_recv_p
integer, dimension(:,:), allocatable recvstatus_srl
integer, dimension(:), allocatable sendrequest_srl
integer, dimension(:,:), allocatable recvstatus_r
integer, dimension(0:3^d &) sizes_p_send_total
integer, dimension(^nd, 0:3) ixr_p_stg_
integer, dimension(-1:1^d &), target type_recv_srl_p1
integer, dimension(:,:), allocatable sendstatus_r
integer, dimension(0:3^d &), target type_recv_p_f
integer, dimension( :^d &), pointer type_recv_srl
integer, dimension(:,:), allocatable sendstatus_p
integer, dimension(-1:1^d &) sizes_srl_recv_total
integer, dimension(:), allocatable recvrequest_p
integer, dimension(:), allocatable recvrequest_srl
integer, dimension(^nd,-1:1^d &) sizes_srl_send_stg
integer, dimension(:,:), allocatable recvstatus_c_p
This module contains definitions of global parameters and variables and some generic functions/subrou...
logical internalboundary
if there is an internal boundary
integer ixghi
Upper index of grid block arrays.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
double precision, dimension(:), allocatable, parameter d
logical ghost_copy
whether copy values instead of interpolation in ghost cells of finer blocks
integer ierrmpi
A global MPI error return code.
integer nghostcells
Number of ghost cells surrounding a grid.
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_convert), pointer phys_to_primitive
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
procedure(sub_convert), pointer phys_to_conserved
character(len=name_len) physics_type
String describing the physics type of the simulation.