MONC
Functions/Subroutines | Variables
simplesetup_mod Module Reference

Functions/Subroutines

type(component_descriptor_type) function, public simplesetup_get_descriptor ()
 
subroutine initialisation_callback (current_state)
 
subroutine allocate_prognostics (current_state)
 
subroutine allocate_prognostic (field, alloc_z, alloc_y, alloc_x, z_grid, y_grid, x_grid)
 
subroutine decompose_grid (current_state)
 
subroutine create_grid (current_state, specific_grid)
 
subroutine read_configuration (current_state)
 

Variables

integer x_size
 
integer y_size
 
integer z_size
 
real(kind=default_precision) zztop
 
real(kind=default_precision) dxx
 
real(kind=default_precision) dyy
 
logical enable_theta =.false.
 

Function/Subroutine Documentation

◆ allocate_prognostic()

subroutine simplesetup_mod::allocate_prognostic ( type(prognostic_field_type), intent(inout)  field,
integer, intent(in)  alloc_z,
integer, intent(in)  alloc_y,
integer, intent(in)  alloc_x,
integer, intent(in)  z_grid,
integer, intent(in)  y_grid,
integer, intent(in)  x_grid 
)
private

Definition at line 96 of file simplesetup.F90.

96  type(prognostic_field_type), intent(inout) :: field
97  integer, intent(in) :: alloc_z, alloc_y, alloc_x, z_grid, y_grid, x_grid
98 
99  field%active=.true.
100  field%grid(z_index) = z_grid
101  field%grid(y_index) = y_grid
102  field%grid(x_index) = x_grid
103  allocate(field%data(alloc_z, alloc_y, alloc_x))
104  field%data=0.0_default_precision
Here is the caller graph for this function:

◆ allocate_prognostics()

subroutine simplesetup_mod::allocate_prognostics ( type(model_state_type), intent(inout)  current_state)
private

Definition at line 46 of file simplesetup.F90.

46  type(model_state_type), intent(inout) :: current_state
47 
48  integer :: alloc_z, alloc_y, alloc_x, i
49 
50  alloc_z=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
51  alloc_y=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
52  alloc_x=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
53 
54 #ifdef U_ACTIVE
55  call allocate_prognostic(current_state%u, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
56  call allocate_prognostic(current_state%zu, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
57  call allocate_prognostic(current_state%su, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
58  call allocate_prognostic(current_state%savu, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
59 #endif
60 #ifdef V_ACTIVE
61  call allocate_prognostic(current_state%v, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
62  call allocate_prognostic(current_state%zv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
63  call allocate_prognostic(current_state%sv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
64  call allocate_prognostic(current_state%savv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
65 #endif
66 #ifdef W_ACTIVE
67  call allocate_prognostic(current_state%w, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
68  call allocate_prognostic(current_state%zw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
69  call allocate_prognostic(current_state%sw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
70  call allocate_prognostic(current_state%savw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
71 #endif
72 
73  if (enable_theta) then
74  call allocate_prognostic(current_state%th, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
75  call allocate_prognostic(current_state%zth, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
76  call allocate_prognostic(current_state%sth, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
77  end if
78  if (current_state%number_q_fields .gt. 0) then
79  allocate(current_state%q(current_state%number_q_fields), &
80  current_state%zq(current_state%number_q_fields), current_state%sq(current_state%number_q_fields))
81  do i=1, current_state%number_q_fields
82  call allocate_prognostic(current_state%q(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
83  call allocate_prognostic(current_state%zq(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
84  call allocate_prognostic(current_state%sq(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
85  end do
86 
87  ! Set standard q indices
88  current_state%water_vapour_mixing_ratio_index=get_q_index(standard_q_names%VAPOUR, 'simplesetup')
89  current_state%liquid_water_mixing_ratio_index=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'simplesetup')
90  end if
91 
92 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ create_grid()

subroutine simplesetup_mod::create_grid ( type(model_state_type), intent(inout)  current_state,
type(global_grid_type), intent(inout)  specific_grid 
)
private

Definition at line 118 of file simplesetup.F90.

118  type(model_state_type), intent(inout) :: current_state
119  type(global_grid_type), intent(inout) :: specific_grid
120 
121  integer, parameter :: kgd_size=20
122  integer :: number_kgd, i, kgd(kgd_size)
123  real(kind=DEFAULT_PRECISION) :: hgd(kgd_size)
124 
125  kgd=-1
126 
127  call options_get_integer_array(current_state%options_database, "kgd", kgd)
128  call options_get_real_array(current_state%options_database, "hgd", hgd)
129 
130  if (kgd(1)==1)then
131  if (hgd(1)/=0.0_default_precision)then
132  call log_log(log_error, "Lowest level is assumed to lie at the surface, check hgd(1)")
133  else
134  kgd(1:kgd_size-1) = kgd(2:)
135  hgd(1:kgd_size-1) = hgd(2:)
136  end if
137  end if
138 
139  do i=1,size(kgd)
140  if (kgd(i) == -1) exit
141  end do
142  number_kgd=i-1
143 
144  if (number_kgd .gt. 0) then
145  allocate(current_state%global_grid%configuration%vertical%kgd(number_kgd), &
146  current_state%global_grid%configuration%vertical%hgd(number_kgd))
147  current_state%global_grid%configuration%vertical%kgd=kgd(1:number_kgd)
148  current_state%global_grid%configuration%vertical%hgd=hgd(1:number_kgd)
149  end if
150 
151  specific_grid%bottom(z_index) = 0
152  specific_grid%bottom(y_index) = 0
153  specific_grid%bottom(x_index) = 0
154 
155  specific_grid%top(z_index) = zztop
156  specific_grid%top(y_index) = dyy * y_size
157  specific_grid%top(x_index) = dxx * x_size
158 
159  specific_grid%resolution(z_index) = zztop / z_size
160  specific_grid%resolution(y_index) = dyy
161  specific_grid%resolution(x_index) = dxx
162 
163  specific_grid%size(z_index) = z_size
164  specific_grid%size(y_index) = y_size
165  specific_grid%size(x_index) = x_size
166 
167  specific_grid%active(z_index) = .true.
168  specific_grid%active(y_index) = .true.
169  specific_grid%active(x_index) = .true.
170 
171  specific_grid%dimensions = 3
Here is the caller graph for this function:

◆ decompose_grid()

subroutine simplesetup_mod::decompose_grid ( type(model_state_type), intent(inout)  current_state)
private

Definition at line 108 of file simplesetup.F90.

108  type(model_state_type), intent(inout) :: current_state
109 
110  if (associated(current_state%parallel%decomposition_procedure)) then
111  call current_state%parallel%decomposition_procedure(current_state)
112  else
113  call log_log(log_error, "No decomposition specified")
114  end if
Here is the caller graph for this function:

◆ initialisation_callback()

subroutine simplesetup_mod::initialisation_callback ( type(model_state_type), intent(inout), target  current_state)
private

Definition at line 30 of file simplesetup.F90.

30  type(model_state_type), intent(inout), target :: current_state
31 
32  call read_configuration(current_state)
33  if (.not. current_state%initialised) then
34  current_state%dtm=options_get_real(current_state%options_database, "dtm")
35  current_state%dtm_new=current_state%dtm
36  call create_grid(current_state, current_state%global_grid)
37  call decompose_grid(current_state)
38  call allocate_prognostics(current_state)
39  current_state%initialised=.true.
40 
41  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_configuration()

subroutine simplesetup_mod::read_configuration ( type(model_state_type), intent(inout), target  current_state)
private

Definition at line 175 of file simplesetup.F90.

175  type(model_state_type), intent(inout), target :: current_state
176 
177  current_state%rhobous=options_get_real(current_state%options_database, "rhobous")
178  current_state%thref0=options_get_real(current_state%options_database, "thref0")
179  current_state%number_q_fields=options_get_integer(current_state%options_database, "number_q_fields")
180  current_state%surface_pressure=options_get_real(current_state%options_database, "surface_pressure")
181  current_state%surface_reference_pressure=options_get_real(current_state%options_database, "surface_reference_pressure")
182 
183  current_state%use_anelastic_equations=options_get_logical(current_state%options_database, "use_anelastic_equations")
184 
185  current_state%origional_vertical_grid_setup=options_get_logical(current_state%options_database, &
186  "origional_vertical_grid_setup")
187  current_state%passive_q=options_get_logical(current_state%options_database, "passive_q")
188  current_state%passive_th=options_get_logical(current_state%options_database, "passive_th")
189  current_state%rmlmax=options_get_real(current_state%options_database, "rmlmax")
190  current_state%calculate_th_and_q_init=options_get_logical(current_state%options_database, "calculate_th_and_q_init")
191  current_state%use_viscosity_and_diffusion=options_get_logical(current_state%options_database, "use_viscosity_and_diffusion")
192  current_state%backscatter=options_get_logical(current_state%options_database, "backscatter")
193 
194  x_size=options_get_integer(current_state%options_database, "x_size")
195  y_size=options_get_integer(current_state%options_database, "y_size")
196  z_size=options_get_integer(current_state%options_database, "z_size")
197  dxx=options_get_real(current_state%options_database, "dxx")
198  dyy=options_get_real(current_state%options_database, "dyy")
199  zztop=options_get_real(current_state%options_database, "zztop")
200  enable_theta=options_get_logical(current_state%options_database, "enable_theta")
201 
202  if (current_state%rmlmax<=0.0)current_state%rmlmax=0.23 * max(dxx, dyy)
203 
204  if (.not. enable_theta) current_state%passive_th=.true.
205  if ( current_state%number_q_fields == 0) current_state%passive_q=.true.
206 
207  current_state%galilean_transformation=options_get_logical(current_state%options_database, "galilean_transformation")
208  if (current_state%galilean_transformation)then
209  current_state%fix_ugal=options_get_logical(current_state%options_database, "fix_ugal")
210  current_state%fix_vgal=options_get_logical(current_state%options_database, "fix_vgal")
211  if (current_state%fix_ugal)current_state%ugal=options_get_real(current_state%options_database, "ugal")
212  if (current_state%fix_vgal)current_state%vgal=options_get_real(current_state%options_database, "vgal")
213  end if
214 
Here is the caller graph for this function:

◆ simplesetup_get_descriptor()

type(component_descriptor_type) function, public simplesetup_mod::simplesetup_get_descriptor ( )

Definition at line 24 of file simplesetup.F90.

24  simplesetup_get_descriptor%name="simplesetup"
25  simplesetup_get_descriptor%version=0.1
26  simplesetup_get_descriptor%initialisation=>initialisation_callback
Here is the call graph for this function:

Variable Documentation

◆ dxx

real(kind=default_precision) simplesetup_mod::dxx
private

Definition at line 18 of file simplesetup.F90.

◆ dyy

real(kind=default_precision) simplesetup_mod::dyy
private

Definition at line 18 of file simplesetup.F90.

◆ enable_theta

logical simplesetup_mod::enable_theta =.false.
private

Definition at line 19 of file simplesetup.F90.

19  logical :: enable_theta=.false.

◆ x_size

integer simplesetup_mod::x_size
private

Definition at line 17 of file simplesetup.F90.

17  integer :: x_size, y_size, z_size

◆ y_size

integer simplesetup_mod::y_size
private

Definition at line 17 of file simplesetup.F90.

◆ z_size

integer simplesetup_mod::z_size
private

Definition at line 17 of file simplesetup.F90.

◆ zztop

real(kind=default_precision) simplesetup_mod::zztop
private

Definition at line 18 of file simplesetup.F90.

18  real(kind=DEFAULT_PRECISION) :: zztop, dxx, dyy