MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_coarsen.t
Go to the documentation of this file.
1 module mod_coarsen
2 
3  implicit none
4  private
5 
6  public :: coarsen_grid
7 
8 contains
9 
10 
11  !> coarsen one grid to its coarser representative
12  subroutine coarsen_grid(sFi,ixFiG^L,ixFi^L,sCo,ixCoG^L,ixCo^L)
14  use mod_physics
15 
16  type(state), intent(inout) :: sfi, sco
17  integer, intent(in) :: ixfig^l, ixfi^l, ixcog^l, ixco^l
18 
19  double precision :: cofiratio
20  double precision :: b_energy_change(ixcog^s)
21  integer :: ixco^d, ixfi^d, iw
22 
23  associate(wfi=>sfi%w(ixfig^s,1:nw), wco=>sco%w(ixcog^s,1:nw))
24  staggered: associate(wfis=>sfi%ws,wcos=>sco%ws)
25  ! coarsen by 2 in every direction - conservatively
26 
27  if(coarsenprimitive) call phys_to_primitive(ixfig^l,ixfi^l,wfi,sfi%x)
28 
29  if(slab_uniform) then
30  cofiratio=one/dble(2**ndim)
31  do iw=1,nw
32  {do ixco^db = ixco^lim^db
33  ixfi^db=2*(ixco^db-ixcomin^db)+ixfimin^db\}
34  wco(ixco^d,iw)=sum(wfi(ixfi^d:ixfi^d+1,iw))*cofiratio
35  {end do\}
36  end do
37  else
38  do iw=1,nw
39  {do ixco^db = ixco^lim^db
40  ixfi^db=2*(ixco^db-ixcomin^db)+ixfimin^db\}
41  wco(ixco^d,iw)= &
42  sum(sfi%dvolume(ixfi^d:ixfi^d+1)*wfi(ixfi^d:ixfi^d+1,iw)) &
43  /sco%dvolume(ixco^d)
44  {end do\}
45  end do
46  end if
47 
48  if(stagger_grid) then
49  do iw=1,nws
50  ! Start one layer before
51  {do ixco^db = ixcomin^db-kr(^db,iw),ixcomax^db
52  ixfi^db=2*(ixco^db-ixcomin^db+kr(^db,iw))+ixfimin^db-kr(^db,iw)\}
53  ! This if statement catches the axis where surface is zero:
54  if (sco%surfaceC(ixco^d,iw)>1.0d-9*sco%dvolume(ixco^d)) then ! Normal case
55  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)) &
56  /sco%surfaceC(ixco^d,iw)
57  else ! On axis
58  wcos(ixco^d,iw)=zero
59  end if
60  {end do\}
61  end do
62  if(phys_total_energy.and. .not.coarsenprimitive) then
63  b_energy_change(ixco^s)=0.5d0*sum(wco(ixco^s,iw_mag(:))**2,dim=ndim+1)
64  end if
65  ! average to fill cell-centred values
66  call phys_face_to_center(ixco^l,sco)
67  if(phys_total_energy.and. .not.coarsenprimitive) then
68  wco(ixco^s,iw_e)=wco(ixco^s,iw_e)-b_energy_change(ixco^s)+&
69  0.5d0*sum(wco(ixco^s,iw_mag(:))**2,dim=ndim+1)
70  end if
71  end if
72 
73  if(coarsenprimitive) then
74  call phys_to_conserved(ixfig^l,ixfi^l,wfi,sfi%x)
75  call phys_to_conserved(ixcog^l,ixco^l,wco,sco%x)
76  end if
77  end associate staggered
78  end associate
79  end subroutine coarsen_grid
80 
81 end module mod_coarsen
subroutine, public coarsen_grid(sFi, ixFiGL, ixFiL, sCo, ixCoGL, ixCoL)
coarsen one grid to its coarser representative
Definition: mod_coarsen.t:13
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 coarsenprimitive
coarsen primitive variables in level-jump ghost cells
double precision, dimension(:), allocatable, parameter d
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