15 integer :: igrid, iigrid, ixcog^
l
16 double precision :: factor
17 logical,
dimension(:,:),
allocatable :: refine2
19 if (igridstail==0)
return
27 do iigrid=1,igridstail; igrid=igrids(iigrid);
35 do iigrid=1,igridstail; igrid=igrids(iigrid);
41 call mpistop(
"Unknown error estimator")
52 do iigrid=1,igridstail; igrid=igrids(iigrid);
69 integer,
intent(in) :: igrid
71 integer :: iflag, idims, idims2, level
72 integer :: ix^L, hx^L, jx^L, h2x^L, j2x^L, ix^D
73 double precision :: epsilon, threshold, wtol(1:nw), xtol(1:ndim)
74 double precision,
dimension(ixM^T) :: numerator, denominator, error
75 double precision,
dimension(ixG^T) :: tmp, tmp1, tmp2
76 double precision :: w(ixG^T,1:nw)
77 logical,
dimension(ixG^T) :: refineflag, coarsenflag
85 w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
89 w(ixg^t,iw_e)=w(ixg^t,iw_e)+0.5d0*sum(ps(igrid)%B0(ixg^t,:,0)**2,dim=ndim+1) &
90 + sum(w(ixg^t,iw_mag(:))*ps(igrid)%B0(ixg^t,:,0),dim=ndim+1)
91 w(ixg^t,iw_mag(:))=w(ixg^t,iw_mag(:))+ps(igrid)%B0(ixg^t,:,0)
101 call mpistop(
"usr_var_for_errest not defined")
108 hx^l=ix^l-
kr(^d,idims);
109 jx^l=ix^l+
kr(^d,idims);
112 tmp(ix^s)=dlog10(w(jx^s,iflag))-dlog10(w(hx^s,iflag))
114 tmp(ix^s)=w(jx^s,iflag)-w(hx^s,iflag)
118 tmp(ix^s)=dlog10(tmp1(jx^s))-dlog10(tmp1(hx^s))
120 tmp(ix^s)=tmp1(jx^s)-tmp1(hx^s)
126 numerator=numerator+(tmp(j2x^s)-tmp(h2x^s))**2
133 tmp=dabs(dlog10(w(ixg^t,iflag)))
135 tmp=dabs(w(ixg^t,iflag))
139 tmp=dabs(dlog10(tmp1(ixg^t)))
141 tmp=dabs(tmp1(ixg^t))
144 hx^l=ix^l-
kr(^d,idims);
145 jx^l=ix^l+
kr(^d,idims);
146 tmp2(ix^s)=tmp(jx^s)+tmp(hx^s)
151 tmp(
ixm^t)=dabs(dlog10(w(jx^s,iflag))&
152 -dlog10(w(
ixm^t,iflag))) &
153 +dabs(dlog10(w(
ixm^t,iflag))&
154 -dlog10(w(hx^s,iflag)))
156 tmp(
ixm^t)=dabs(w(jx^s,iflag)-w(
ixm^t,iflag)) &
157 +dabs(w(
ixm^t,iflag)-w(hx^s,iflag))
161 tmp(
ixm^t)=dabs(dlog10(tmp1(jx^s))-dlog10(tmp1(
ixm^t))) &
162 +dabs(dlog10(tmp1(
ixm^t))-dlog10(tmp1(hx^s)))
164 tmp(
ixm^t)=dabs(tmp1(jx^s)-tmp1(
ixm^t)) &
165 +dabs(tmp1(
ixm^t)-tmp1(hx^s))
171 denominator=denominator &
175 error=error+
w_refine_weight(iflag)*dsqrt(numerator/max(denominator,epsilon))
181 {
do ix^db=ixmlo^db,ixmhi^db\}
184 wtol(1:nw) = w(ix^d,1:nw)
185 xtol(1:ndim) = ps(igrid)%x(ix^d,1:ndim)
189 if (error(ix^d) >= threshold)
then
190 refineflag(ix^d) = .true.
192 coarsenflag(ix^d) = .true.
196 if (any(refineflag(ixm^t)).and.level<refine_max_level) refine(igrid,mype)=.true.
197 if (all(coarsenflag(ixm^t)).and.level>1) coarsen(igrid,mype)=.true.
206 integer,
intent(in) :: igrid
208 integer :: iflag, idims, level
209 integer :: ix^L, hx^L, jx^L, ix^D
210 double precision :: epsilon, threshold, wtol(1:nw), xtol(1:ndim)
211 double precision,
dimension(ixM^T) :: numerator, denominator, error
212 double precision,
dimension(ixG^T) :: dp, dm, dref, tmp1
213 logical,
dimension(ixG^T) :: refineflag, coarsenflag
227 call mpistop(
"usr_var_for_errest not defined")
234 hx^l=ix^l-
kr(^d,idims);
235 jx^l=ix^l+
kr(^d,idims);
238 dp(ix^s)=dlog10(ps(igrid)%w(jx^s,iflag))-dlog10(ps(igrid)%w(ix^s,iflag))
239 dm(ix^s)=dlog10(ps(igrid)%w(ix^s,iflag))-dlog10(ps(igrid)%w(hx^s,iflag))
240 dref(
ixm^t)=dabs(dlog10(ps(igrid)%w(jx^s,iflag)))&
241 + 2.0d0 * dabs(dlog10(ps(igrid)%w(
ixm^t,iflag))) &
242 + dabs(dlog10(ps(igrid)%w(hx^s,iflag)))
244 dp(ix^s)=ps(igrid)%w(jx^s,iflag)-ps(igrid)%w(ix^s,iflag)
245 dm(ix^s)=ps(igrid)%w(ix^s,iflag)-ps(igrid)%w(hx^s,iflag)
246 dref(
ixm^t)=dabs(ps(igrid)%w(jx^s,iflag))+2.0d0*dabs(ps(igrid)%w(
ixm^t,iflag)) &
247 +dabs(ps(igrid)%w(hx^s,iflag))
251 dp(ix^s)=dlog10(tmp1(jx^s))-dlog10(tmp1(ix^s))
252 dm(ix^s)=dlog10(tmp1(ix^s))-dlog10(tmp1(hx^s))
253 dref(ix^s)=dabs(dlog10(tmp1(jx^s)))&
254 + 2.0d0 * dabs(dlog10(tmp1(ix^s))) &
255 + dabs(dlog10(tmp1(hx^s)))
257 dp(ix^s)=tmp1(jx^s)-tmp1(ix^s)
258 dm(ix^s)=tmp1(ix^s)-tmp1(hx^s)
259 dref(ix^s)=dabs(tmp1(jx^s))+2.0d0*dabs(tmp1(ix^s)) &
264 numerator(
ixm^t)=numerator+(dp(
ixm^t)-dm(
ixm^t))**2
265 denominator(
ixm^t)=denominator &
269 error=error+
w_refine_weight(iflag)*dsqrt(numerator/max(denominator,epsilon))
276 {
do ix^db=ixmlo^db,ixmhi^db\}
279 wtol(1:nw) = ps(igrid)%w(ix^d,1:nw)
280 xtol(1:ndim) = ps(igrid)%x(ix^d,1:ndim)
284 if (error(ix^d) >= threshold)
then
285 refineflag(ix^d) = .true.
287 coarsenflag(ix^d) = .true.
291 if (any(refineflag(ixm^t)).and.level<refine_max_level) refine(igrid,mype)=.true.
292 if (all(coarsenflag(ixm^t)).and.level>1) coarsen(igrid,mype)=.true.
301 integer,
intent(in) :: igrid
302 double precision,
intent(in) :: w(ixG^T,nw)
305 integer :: my_refine, my_coarsen
306 double precision :: qt
307 logical,
dimension(ixG^T) :: refineflag
323 my_refine,my_coarsen)
326 if (my_coarsen==1)
then
336 if (my_coarsen==-1)
then
340 if (my_refine==1)
then
350 if (my_refine==-1)
then
356 refineflag(
ixm^t)=.true.
367 integer,
intent(in) :: igrid
368 double precision,
intent(in) :: w(ixg^t,nw)
370 integer :: level, my_levmin, my_levmax
371 logical,
dimension(ixG^T) :: refineflag
383 if (level>my_levmax)
then
386 elseif (level<my_levmin)
then
391 if (level==my_levmin .or. level==my_levmax)
then
405 integer,
intent(in) :: igrid
406 logical,
dimension(ixG^T),
intent(in) :: refineflag
408 integer :: ishiftbuf^D, i^D, ix^L, ineighbor, ipe_neighbor, level
410 ishiftbuf^d=ixmhi^d-ixmlo^d-
nbufferx^d+1;
412 ixmin^d=max(ixmlo^d,ixmlo^d+i^d*ishiftbuf^d);
413 ixmax^d=min(ixmhi^d,ixmhi^d+i^d*ishiftbuf^d);
414 if (ixmax^d<ixmin^d|.or.) cycle
415 if (any(refineflag(ix^s)))
then
416 select case (neighbor_type(i^d,igrid))
417 case (neighbor_coarse)
418 ineighbor=neighbor(1,i^d,igrid)
419 ipe_neighbor=neighbor(2,i^d,igrid)
420 if (.not.
refine(ineighbor,ipe_neighbor))
then
421 buffer(ineighbor,ipe_neighbor)=.true.
422 refine(ineighbor,ipe_neighbor)=.true.
424 case (neighbor_sibling)
427 ineighbor=neighbor(1,i^d,igrid)
428 ipe_neighbor=neighbor(2,i^d,igrid)
429 if (.not.
refine(ineighbor,ipe_neighbor))
then
430 buffer(ineighbor,ipe_neighbor)=.true.
431 refine(ineighbor,ipe_neighbor)=.true.
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
subroutine refinebuffer(igrid, refineflag)
subroutine forcedrefine_grid(igrid, w)
subroutine lohner_grid(igrid)
subroutine lohner_orig_grid(igrid)
subroutine, public errest
Do all local error estimation which determines (de)refinement.
subroutine, public forcedrefine_grid_io(igrid, w)
Module with basic grid data structures.
logical, dimension(:,:), allocatable, save refine
logical, dimension(:,:), allocatable, save buffer
logical, dimension(:,:), allocatable, save coarsen
AMR flags and grids-in-use identifier per processor (igrid,ipe)
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
double precision global_time
The global simulation time.
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical, dimension(:), allocatable logflag
double precision, dimension(:), allocatable amr_wavefilter
refinement: lohner estimate wavefilter setting
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, parameter plevel_
integer, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
logical b0field
split magnetic field as background B0 field
integer nbufferx
Number of cells as buffer zone.
double precision, dimension(:), allocatable w_refine_weight
Weights of variables used to calculate error for mesh refinement.
double precision, dimension(:), allocatable refine_threshold
Error tolerance for refinement decision.
integer refine_max_level
Maximal number of AMR levels.
double precision, dimension(:), allocatable derefine_ratio
integer max_blocks
The maximum number of grid blocks in a processor.
integer, dimension(:,:), allocatable node
This module defines the procedures of a physics module. It contains function pointers for the various...
logical phys_energy
Solve energy equation or not.
Module with all the methods that users can customize in AMRVAC.
procedure(a_refine_threshold), pointer usr_refine_threshold
procedure(refine_grid), pointer usr_refine_grid
procedure(var_for_errest), pointer usr_var_for_errest