MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_physicaldata.t
Go to the documentation of this file.
2  implicit none
3  save
4 
5  type state
6  !> ID of a grid block
7  integer :: igrid=-1
8  !> index range of block array in cell centers
9  integer :: ixg^l
10  !> index range of block array in cell faces
11  integer :: ixgs^l
12  !> level of AMR
13  integer :: level
14  !> If it face a physical boundary
15  logical, dimension(:), pointer :: is_physical_boundary(:) =>null()
16  !> Variables, normally cell center conservative values
17  double precision, dimension(:^D&,:), allocatable :: w
18  !> Variables, cell face values
19  double precision, dimension(:^D&,:), allocatable :: ws
20  !> Variables, cell edge values
21  double precision, dimension(:^D&,:), allocatable :: we
22  !> Variables, cell corner values
23  double precision, dimension(:^D&,:), allocatable :: wc
24  !> extra variables do not need ghost cell and equation flux
25  double precision, dimension(:^D&,:), pointer :: wextra=>null()
26  !> Time-independent magnetic field at cell center and cell interface
27  double precision, dimension(:^D&,:,:), pointer :: b0=>null()
28  !> Time-independent electric current density at cell center
29  double precision, dimension(:^D&,:), pointer :: j0=>null()
30  !> Time-independent equi vars (B0 is not into this array)
31  double precision, dimension(:^D&,:,:), pointer :: equi_vars=>null()
32  !> Cell-center positions
33  double precision, dimension(:^D&,:), pointer :: x=>null()
34  !> Cell sizes in coordinate units
35  double precision, dimension(:^D&,:), pointer :: dx=>null()
36  !> Cell local timesteps
37  double precision, dimension(:^D&), pointer :: dt=>null()
38  !> Cell sizes at cell center in length unit
39  double precision, dimension(:^D&,:), pointer :: ds=>null()
40  !> Cell sizes at cell face in length unit
41  double precision, dimension(:^D&,:), pointer :: dsc=>null()
42  !> Volumes of a cell
43  double precision, dimension(:^D&), pointer :: dvolume=>null()
44  !> Areas of cell-center surfaces
45  double precision, dimension(:^D&,:), pointer :: surface=>null()
46  !> Areas of cell-face surfaces
47  double precision, dimension(:^D&,:), pointer :: surfacec=>null()
48  !> special values for a block
49  double precision, dimension(:), pointer :: special_values=>null()
50  end type state
51 
52 {^nooned
53  type state_sub
54  !> ID of a grid block
55  integer :: igrid=-1
56  !> Variables, normally center
57  double precision, dimension(:^DE&,:), allocatable :: w
58  !> Variables for the cornerpositions on the slice
59  double precision, dimension(:^DE&,:), allocatable :: wc
60  !> Variables, normally center, one level coarser representative
61  double precision, dimension(:^DE&,:), allocatable :: wcoarse
62  !> Cell-center positions
63  double precision, dimension(:^DE&,:), allocatable :: x
64  !> Corner positions on the slice
65  double precision, dimension(:^DE&,:), allocatable :: xc
66  !> Cell-center positions, one level coarser representative
67  double precision, dimension(:^DE&,:), allocatable :: xcoarse
68  !> Cell sizes
69  double precision, dimension(:^DE&,:), allocatable :: dx
70  !> Cell sizes, one level coarser
71  double precision, dimension(:^D&,:), allocatable :: dxcoarse
72  !> Cell sizes in length unit
73  double precision, dimension(:^D&,:), allocatable :: ds
74  !> Volumes of a cell
75  double precision, dimension(:^DE&), allocatable :: dvolume
76  !> Volumes of a cell, one level coarser representative
77  double precision, dimension(:^DE&), allocatable :: dvolumecoarse
78  !> Areas of cell-center surfaces
79  double precision, dimension(:^DE&,:), allocatable :: surface
80  !> Areas of cell-face surfaces
81  double precision, dimension(:^DE&,:), allocatable :: surfacec
82  end type state_sub
83 }
84 {^ifoned
85  type state_sub
86  !> ID of a grid block
87  integer :: igrid=-1
88  !> Variables, normally center
89  double precision, dimension(:), allocatable :: w
90  !> Variables for the cornerpositions on the slice
91  double precision, dimension(:), allocatable :: wc
92  !> Variables, normally center, one level coarser representative
93  double precision, dimension(:), allocatable :: wcoarse
94  !> Cell-center positions
95  double precision, dimension(:), allocatable :: x
96  !> Corner positions on the slice
97  double precision, dimension(:), allocatable :: xc
98  !> Cell-center positions, one level coarser representative
99  double precision, dimension(:), allocatable :: xcoarse
100  end type state_sub
101 }
103  !> Variables new state
104  double precision, dimension(:^D&,:), allocatable :: w
105  !> Variables old state
106  double precision, dimension(:^D&,:), allocatable :: wold
107  end type grid_field
108  !> buffer for pole boundary
109  type(state) :: pole_buf
110 
111  !> array of physical states for all blocks on my processor
112  type(state), dimension(:), allocatable, target :: ps
113  !> array of physical states, temp 1 for multi-step time integrator
114  type(state), dimension(:), allocatable, target :: ps1
115  !> array of physical states, temp 2 for multi-step time integrator
116  type(state), dimension(:), allocatable, target :: ps2
117  !> array of physical states, temp 3 for multi-step time integrator
118  type(state), dimension(:), allocatable, target :: ps3
119  !> array of physical states, temp 4 for multi-step time integrator
120  type(state), dimension(:), allocatable, target :: ps4
121  !> array of physical blocks, one level coarser representative
122  type(state), dimension(:), allocatable, target :: psc
123 
124  !> array of physical blocks in reduced dimension
125  type(state_sub), dimension(:), allocatable, target :: ps_sub
126 
127 {^ifoned
128  double precision, dimension(:), allocatable :: collapseddata
129 }
130 {^nooned
131  double precision, dimension(:^DE&,:), allocatable :: collapseddata
132 }
133  !> array of physical blocks of meshed fields for particles
134  type(grid_field), dimension(:), allocatable, target :: gridvars
135 
136  !> velocities store for constrained transport
138  double precision, dimension(:^D&,:), allocatable :: vnorm,cbarmin,cbarmax
139  double precision, dimension(:^D&,:,:), allocatable :: vbarc,vbarlc,vbarrc
140  end type ct_velocity
141 
142 end module mod_physicaldata
type(state), dimension(:), allocatable, target ps4
array of physical states, temp 4 for multi-step time integrator
type(state_sub), dimension(:), allocatable, target ps_sub
array of physical blocks in reduced dimension
type(grid_field), dimension(:), allocatable, target gridvars
array of physical blocks of meshed fields for particles
type(state), dimension(:), allocatable, target ps2
array of physical states, temp 2 for multi-step time integrator
type(state), dimension(:), allocatable, target ps
array of physical states for all blocks on my processor
type(state), dimension(:), allocatable, target psc
array of physical blocks, one level coarser representative
type(state), dimension(:), allocatable, target ps1
array of physical states, temp 1 for multi-step time integrator
type(state) pole_buf
buffer for pole boundary
double precision, dimension(:^de &,:), allocatable collapseddata
type(state), dimension(:), allocatable, target ps3
array of physical states, temp 3 for multi-step time integrator
velocities store for constrained transport