2D AMRCLAW
flag2refine2.f90
Go to the documentation of this file.
1 ! ::::::::::::::::::::: flag2refine ::::::::::::::::::::::::::::::::::
2 !
51 !
52 ! ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
53 subroutine flag2refine2(mx,my,mbc,mbuff,meqn,maux,xlower,ylower,dx,dy,t,level, &
54  tolsp,q,aux,amrflags,DONTFLAG,DOFLAG)
55 
56  use regions_module
57 
58  implicit none
59 
60  ! Subroutine arguments
61  integer, intent(in) :: mx,my,mbc,meqn,maux,level,mbuff
62  real(kind=8), intent(in) :: xlower,ylower,dx,dy,t,tolsp
63 
64  real(kind=8), intent(in) :: q(meqn,1-mbc:mx+mbc,1-mbc:my+mbc)
65  real(kind=8), intent(in) :: aux(maux,1-mbc:mx+mbc,1-mbc:my+mbc)
66 
67  ! Flagging
68  real(kind=8),intent(inout) :: amrflags(1-mbuff:mx+mbuff,1-mbuff:my+mbuff)
69  real(kind=8), intent(in) :: DONTFLAG
70  real(kind=8), intent(in) :: DOFLAG
71 
72  logical :: allowflag
73  external allowflag
74 
75  ! Locals
76  integer :: i,j,m
77  real(kind=8) :: x_c,y_c,x_low,y_low,x_hi,y_hi
78  real(kind=8) :: dqi(meqn), dqj(meqn), dq(meqn)
79 
80  ! Initialize flags
81  amrflags = dontflag
82 
83  ! Loop over interior points on this grid
84  ! (i,j) grid cell is [x_low,x_hi] x [y_low,y_hi], cell center at (x_c,y_c)
85  ! This information is not needed for the default flagging based on
86  ! undivided differences, but might be needed in a user's version.
87  ! Note that if you want to refine only in certain space-time regions,
88  ! it may be easiest to use the "regions" feature. The flags set here or
89  ! in the Richardson error estimator are potentially modified by the
90  ! min_level and max_level specified in any regions.
91 
92  y_loop: do j=1,my
93  y_low = ylower + (j - 1) * dy
94  y_c = ylower + (j - 0.5d0) * dy
95  y_hi = ylower + j * dy
96 
97  x_loop: do i = 1,mx
98  x_low = xlower + (i - 1) * dx
99  x_c = xlower + (i - 0.5d0) * dx
100  x_hi = xlower + i * dx
101 
102  ! -----------------------------------------------------------------
103  dq = 0.d0
104  dqi = abs(q(:,i+1,j) - q(:,i-1,j))
105  dqj = abs(q(:,i,j+1) - q(:,i,j-1))
106  dq = max(dq,dqi,dqj)
107 
108  ! default checks all components of undivided difference:
109  do m=1,meqn
110  if (dq(m) > tolsp) then
111  amrflags(i,j) = doflag
112  cycle x_loop
113  endif
114  enddo
115 
116  enddo x_loop
117  enddo y_loop
118 
119 end subroutine flag2refine2
subroutine flag2refine2(mx, my, mbc, mbuff, meqn, maux, xlower, ylower, dx, dy, t, level, tolsp, q, aux, amrflags, DONTFLAG, DOFLAG)
User routine to control flagging of points for refinement.