MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
coarsen.t
Go to the documentation of this file.
1 !> coarsen one grid to its coarser representative
2 subroutine coarsen_grid(sFi,ixFiG^L,ixFi^L,sCo,ixCoG^L,ixCo^L)
4  use mod_physics
5 
6  type(state), intent(inout) :: sFi, sCo
7  integer, intent(in) :: ixFiG^L, ixFi^L, ixCoG^L, ixCo^L
8 
9  integer :: ixCo^D, ixFi^D, iw
10  double precision :: CoFiratio
11  double precision :: B_energy_change(ixCoG^S)
12 
13  associate(wfi=>sfi%w(ixfig^s,1:nw), wco=>sco%w(ixcog^s,1:nw))
14  staggered: associate(wfis=>sfi%ws,wcos=>sco%ws)
15  ! coarsen by 2 in every direction - conservatively
16 
17  if(coarsenprimitive) call phys_to_primitive(ixfig^l,ixfi^l,wfi,sfi%x)
18 
19  if(slab_uniform) then
20  cofiratio=one/dble(2**ndim)
21  do iw=1,nw
22  {do ixco^db = ixco^lim^db
23  ixfi^db=2*(ixco^db-ixcomin^db)+ixfimin^db\}
24  wco(ixco^d,iw)=sum(wfi(ixfi^d:ixfi^d+1,iw))*cofiratio
25  {end do\}
26  end do
27  else
28  do iw=1,nw
29  {do ixco^db = ixco^lim^db
30  ixfi^db=2*(ixco^db-ixcomin^db)+ixfimin^db\}
31  wco(ixco^d,iw)= &
32  sum(sfi%dvolume(ixfi^d:ixfi^d+1)*wfi(ixfi^d:ixfi^d+1,iw)) &
33  /sco%dvolume(ixco^d)
34  {end do\}
35  end do
36  end if
37 
38  if(stagger_grid) then
39  do iw=1,nws
40  ! Start one layer before
41  {do ixco^db = ixcomin^db-kr(^db,iw),ixcomax^db
42  ixfi^db=2*(ixco^db-ixcomin^db+kr(^db,iw))+ixfimin^db-kr(^db,iw)\}
43  ! This if statement catches the axis where surface is zero:
44  if (sco%surfaceC(ixco^d,iw)>1.0d-9*sco%dvolume(ixco^d)) then ! Normal case
45  wcos(ixco^d,iw)=sum(sfi%surfaceC(ixfi^d:ixfi^d+1-kr(iw,^d),iw)*wfis(ixfi^d:ixfi^d+1-kr(iw,^d),iw)) &
46  /sco%surfaceC(ixco^d,iw)
47  else ! On axis
48  wcos(ixco^d,iw)=zero
49  end if
50  {end do\}
51  end do
52  if(phys_total_energy.and. .not.coarsenprimitive) then
53  b_energy_change(ixco^s)=0.5d0*sum(wco(ixco^s,iw_mag(:))**2,dim=ndim+1)
54  end if
55  ! average to fill cell-centred values
56  call phys_face_to_center(ixco^l,sco)
57  if(phys_total_energy.and. .not.coarsenprimitive) then
58  wco(ixco^s,iw_e)=wco(ixco^s,iw_e)-b_energy_change(ixco^s)+&
59  0.5d0*sum(wco(ixco^s,iw_mag(:))**2,dim=ndim+1)
60  end if
61  end if
62 
63  if(coarsenprimitive) then
64  call phys_to_conserved(ixfig^l,ixfi^l,wfi,sfi%x)
65  call phys_to_conserved(ixcog^l,ixco^l,wco,sco%x)
66  end if
67  end associate staggered
68  end associate
69 end subroutine coarsen_grid
subroutine coarsen_grid(sFi, ixFiGL, ixFiL, sCo, ixCoGL, ixCoL)
coarsen one grid to its coarser representative
Definition: coarsen.t:3
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59