MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_initialize_amr.t
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: initlevelone
7  public :: modify_ic
8  public :: initial_condition
10 
11 
12 contains
13 
14 
15  !> Generate and initialize all grids at the coarsest level (level one)
16  subroutine initlevelone
22 
23  integer :: iigrid, igrid{#ifdef evolvingboundary , morton_no}
24 
25  levmin=1
26  levmax=1
27 
28  call init_forest_root
29 
30  call getigrids
32 
33  ! fill solution space of all root grids
34  do iigrid=1,igridstail; igrid=igrids(iigrid);
35  call alloc_node(igrid)
36  ! in case gradient routine used in initial condition, ensure geometry known
37  call initial_condition(igrid)
38  end do
39  {#IFDEF EVOLVINGBOUNDARY
40  ! mark physical-boundary blocks on space-filling curve
41  do morton_no=morton_start(mype),morton_stop(mype)
42  igrid=sfc_to_igrid(morton_no)
43  if (phyboundblock(igrid)) sfc_phybound(morton_no)=1
44  end do
45  call mpi_allreduce(mpi_in_place,sfc_phybound,nleafs,mpi_integer,&
46  mpi_sum,icomm,ierrmpi)
47  }
48 
49  ! update ghost cells
50  call getbc(global_time,0.d0,ps,iwstart,nwgc)
51 
52  end subroutine initlevelone
53 
54  !> fill in initial condition
55  subroutine initial_condition(igrid)
56  ! Need only to set the mesh values (can leave ghost cells untouched)
59  use mod_comm_lib, only: mpistop
60 
61  integer, intent(in) :: igrid
62 
63  ! in case gradient routine used in initial condition, ensure geometry known
64  block=>ps(igrid)
65  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
66 
67  if (.not. associated(usr_init_one_grid)) then
68  call mpistop("usr_init_one_grid not defined")
69  else
70  call usr_init_one_grid(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x)
71  end if
72 
73  end subroutine initial_condition
74 
75  !> modify initial condition
76  subroutine modify_ic
79  use mod_comm_lib, only: mpistop
80 
81  integer :: iigrid, igrid
82 
83  do iigrid=1,igridstail; igrid=igrids(iigrid);
84  block=>ps(igrid)
85  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
86 
87  if (.not. associated(usr_init_one_grid)) then
88  call mpistop("usr_init_one_grid not defined")
89  else
90  call usr_init_one_grid(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x)
91  end if
92  end do
93 
94  end subroutine modify_ic
95 
96  !> improve initial condition after initialization
99  use mod_usr_methods
102  use mod_physics
104 
105  logical :: active
106 
107  if(associated(usr_improve_initial_condition)) then
109  else if(stagger_grid) then
110  if(associated(usr_init_vector_potential)) then
111  ! re-calculate magnetic field from the vector potential in a
112  ! completely divergency free way for AMR mesh in 3D
113  if(levmax>levmin.and.ndim==3) call recalculateb
114  end if
115  if(slab_uniform.and.associated(phys_clean_divb)) then
116  ! Project out the divB using multigrid poisson solver
117  ! if not initialised from vector potential
118  if(.not.use_multigrid) call mg_setup_multigrid()
119  call phys_clean_divb(global_time,0.d0,active)
120  call getbc(global_time,0.d0,ps,iwstart,nwgc)
121  end if
122  end if
123 
124  end subroutine improve_initial_condition
125 
126 end module mod_initialize_amr
subroutine, public alloc_node(igrid)
allocate arrays on igrid node
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
subroutine recalculateb
re-calculate the magnetic field from the vector potential in a completely divergency free way
subroutine, public init_forest_root
build root forest
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
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, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
logical, dimension(:), allocatable phyboundblock
True if a block has any physical boundary.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
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.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
logical use_multigrid
Use multigrid (only available in 2D and 3D)
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision, dimension(ndim) dxlevel
subroutine, public improve_initial_condition()
improve initial condition after initialization
subroutine, public initlevelone
Generate and initialize all grids at the coarsest level (level one)
subroutine, public modify_ic
modify initial condition
subroutine, public initial_condition(igrid)
fill in initial condition
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
subroutine mg_setup_multigrid()
Setup multigrid for usage.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_clean_divb), pointer phys_clean_divb
Definition: mod_physics.t:88
Module with all the methods that users can customize in AMRVAC.
procedure(p_no_args), pointer usr_improve_initial_condition
procedure(init_one_grid), pointer usr_init_one_grid
Initialize earch grid block data.
procedure(init_vector_potential), pointer usr_init_vector_potential