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 100 of file simplesetup.F90.

101  type(prognostic_field_type), intent(inout) :: field
102  integer, intent(in) :: alloc_z, alloc_y, alloc_x, z_grid, y_grid, x_grid
103 
104  field%active=.true.
105  field%grid(z_index) = z_grid
106  field%grid(y_index) = y_grid
107  field%grid(x_index) = x_grid
108  allocate(field%data(alloc_z, alloc_y, alloc_x))
109  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 47 of file simplesetup.F90.

48  type(model_state_type), intent(inout) :: current_state
49 
50  integer :: alloc_z, alloc_y, alloc_x, i
51 
52  alloc_z=current_state%local_grid%size(z_index) + current_state%local_grid%halo_size(z_index) * 2
53  alloc_y=current_state%local_grid%size(y_index) + current_state%local_grid%halo_size(y_index) * 2
54  alloc_x=current_state%local_grid%size(x_index) + current_state%local_grid%halo_size(x_index) * 2
55 
56 #ifdef U_ACTIVE
57  call allocate_prognostic(current_state%u, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
58  call allocate_prognostic(current_state%zu, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
59  call allocate_prognostic(current_state%su, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
60  call allocate_prognostic(current_state%savu, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, primal_grid)
61 #endif
62 #ifdef V_ACTIVE
63  call allocate_prognostic(current_state%v, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
64  call allocate_prognostic(current_state%zv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
65  call allocate_prognostic(current_state%sv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
66  call allocate_prognostic(current_state%savv, alloc_z, alloc_y, alloc_x, dual_grid, primal_grid, dual_grid)
67 #endif
68 #ifdef W_ACTIVE
69  call allocate_prognostic(current_state%w, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
70  call allocate_prognostic(current_state%zw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
71  call allocate_prognostic(current_state%sw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
72  call allocate_prognostic(current_state%savw, alloc_z, alloc_y, alloc_x, primal_grid, dual_grid, dual_grid)
73 #endif
74 
75  if (enable_theta) then
76  call allocate_prognostic(current_state%th, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
77  call allocate_prognostic(current_state%zth, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
78  call allocate_prognostic(current_state%sth, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
79  end if
80  if (current_state%number_q_fields .gt. 0) then
81  allocate(current_state%q(current_state%number_q_fields), &
82  current_state%zq(current_state%number_q_fields), current_state%sq(current_state%number_q_fields))
83  do i=1, current_state%number_q_fields
84  call allocate_prognostic(current_state%q(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
85  call allocate_prognostic(current_state%zq(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
86  call allocate_prognostic(current_state%sq(i), alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
87  end do
88 
89  ! Set standard q indices
90  current_state%water_vapour_mixing_ratio_index=get_q_index(standard_q_names%VAPOUR, 'simplesetup')
91  current_state%liquid_water_mixing_ratio_index=get_q_index(standard_q_names%CLOUD_LIQUID_MASS, 'simplesetup')
92  end if
93 
94  ! Set arrays for radiative heating rates - Note: this should be protected by a switch
95  call allocate_prognostic(current_state%sth_lw, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
96  call allocate_prognostic(current_state%sth_sw, alloc_z, alloc_y, alloc_x, dual_grid, dual_grid, dual_grid)
97 
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 122 of file simplesetup.F90.

123  type(model_state_type), intent(inout) :: current_state
124  type(global_grid_type), intent(inout) :: specific_grid
125 
126  integer, parameter :: KGD_SIZE=200
127  integer :: number_kgd, i, kgd(KGD_SIZE)
128  real(kind=default_precision) :: hgd(kgd_size)
129 
130  kgd=-1
131 
132  call options_get_integer_array(current_state%options_database, "kgd", kgd)
133  call options_get_real_array(current_state%options_database, "hgd", hgd)
134 
135  if (kgd(1)==1)then
136  if (hgd(1)/=0.0_default_precision)then
137  call log_log(log_error, "Lowest level is assumed to lie at the surface, check hgd(1)")
138  else
139  kgd(1:kgd_size-1) = kgd(2:)
140  hgd(1:kgd_size-1) = hgd(2:)
141  end if
142  end if
143 
144  do i=1,size(kgd)
145  if (kgd(i) == -1) exit
146  end do
147  number_kgd=i-1
148 
149  if (number_kgd .gt. 0) then
150  allocate(current_state%global_grid%configuration%vertical%kgd(number_kgd), &
151  current_state%global_grid%configuration%vertical%hgd(number_kgd))
152  current_state%global_grid%configuration%vertical%kgd=kgd(1:number_kgd)
153  current_state%global_grid%configuration%vertical%hgd=hgd(1:number_kgd)
154  end if
155 
156  specific_grid%bottom(z_index) = 0
157  specific_grid%bottom(y_index) = 0
158  specific_grid%bottom(x_index) = 0
159 
160  specific_grid%top(z_index) = zztop
161  specific_grid%top(y_index) = dyy * y_size
162  specific_grid%top(x_index) = dxx * x_size
163 
164  specific_grid%resolution(z_index) = zztop / z_size
165  specific_grid%resolution(y_index) = dyy
166  specific_grid%resolution(x_index) = dxx
167 
168  specific_grid%size(z_index) = z_size
169  specific_grid%size(y_index) = y_size
170  specific_grid%size(x_index) = x_size
171 
172  specific_grid%active(z_index) = .true.
173  specific_grid%active(y_index) = .true.
174  specific_grid%active(x_index) = .true.
175 
176  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 112 of file simplesetup.F90.

113  type(model_state_type), intent(inout) :: current_state
114 
115  if (associated(current_state%parallel%decomposition_procedure)) then
116  call current_state%parallel%decomposition_procedure(current_state)
117  else
118  call log_log(log_error, "No decomposition specified")
119  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 31 of file simplesetup.F90.

32  type(model_state_type), intent(inout), target :: current_state
33 
34  call read_configuration(current_state)
35  if (.not. current_state%initialised) then
36  current_state%dtm=options_get_real(current_state%options_database, "dtm")
37  current_state%dtm_new=current_state%dtm
38  call create_grid(current_state, current_state%global_grid)
39  call decompose_grid(current_state)
40  call allocate_prognostics(current_state)
41  current_state%initialised=.true.
42 
43  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 179 of file simplesetup.F90.

180  type(model_state_type), intent(inout), target :: current_state
181 
182  current_state%rhobous=options_get_real(current_state%options_database, "rhobous")
183  current_state%thref0=options_get_real(current_state%options_database, "thref0")
184  current_state%number_q_fields=options_get_integer(current_state%options_database, "number_q_fields")
185  current_state%surface_pressure=options_get_real(current_state%options_database, "surface_pressure")
186  current_state%surface_reference_pressure=options_get_real(current_state%options_database, "surface_reference_pressure")
187 
188  current_state%use_anelastic_equations=options_get_logical(current_state%options_database, "use_anelastic_equations")
189 
190  current_state%origional_vertical_grid_setup=options_get_logical(current_state%options_database, &
191  "origional_vertical_grid_setup")
192  current_state%passive_q=options_get_logical(current_state%options_database, "passive_q")
193  current_state%passive_th=options_get_logical(current_state%options_database, "passive_th")
194  current_state%rmlmax=options_get_real(current_state%options_database, "rmlmax")
195  current_state%calculate_th_and_q_init=options_get_logical(current_state%options_database, "calculate_th_and_q_init")
196  current_state%use_viscosity_and_diffusion=options_get_logical(current_state%options_database, "use_viscosity_and_diffusion")
197  current_state%backscatter=options_get_logical(current_state%options_database, "backscatter")
198 
199  x_size=options_get_integer(current_state%options_database, "x_size")
200  y_size=options_get_integer(current_state%options_database, "y_size")
201  z_size=options_get_integer(current_state%options_database, "z_size")
202  dxx=options_get_real(current_state%options_database, "dxx")
203  dyy=options_get_real(current_state%options_database, "dyy")
204  zztop=options_get_real(current_state%options_database, "zztop")
205  enable_theta=options_get_logical(current_state%options_database, "enable_theta")
206 
207  if (current_state%rmlmax<=0.0)current_state%rmlmax=0.23 * max(dxx, dyy)
208 
209  if (.not. enable_theta) current_state%passive_th=.true.
210  if ( current_state%number_q_fields == 0) current_state%passive_q=.true.
211 
212  current_state%galilean_transformation=options_get_logical(current_state%options_database, "galilean_transformation")
213  if (current_state%galilean_transformation)then
214  current_state%fix_ugal=options_get_logical(current_state%options_database, "fix_ugal")
215  current_state%fix_vgal=options_get_logical(current_state%options_database, "fix_vgal")
216  if (current_state%fix_ugal)current_state%ugal=options_get_real(current_state%options_database, "ugal")
217  if (current_state%fix_vgal)current_state%vgal=options_get_real(current_state%options_database, "vgal")
218  end if
219 
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 25 of file simplesetup.F90.

26  simplesetup_get_descriptor%name="simplesetup"
27  simplesetup_get_descriptor%version=0.1
28  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 20 of file simplesetup.F90.

◆ dyy

real(kind=default_precision) simplesetup_mod::dyy
private

Definition at line 20 of file simplesetup.F90.

◆ enable_theta

logical simplesetup_mod::enable_theta =.false.
private

Definition at line 21 of file simplesetup.F90.

21  logical :: enable_theta=.false.

◆ x_size

integer simplesetup_mod::x_size
private

Definition at line 19 of file simplesetup.F90.

19  integer :: x_size, y_size, z_size

◆ y_size

integer simplesetup_mod::y_size
private

Definition at line 19 of file simplesetup.F90.

◆ z_size

integer simplesetup_mod::z_size
private

Definition at line 19 of file simplesetup.F90.

◆ zztop

real(kind=default_precision) simplesetup_mod::zztop
private

Definition at line 20 of file simplesetup.F90.

20  real(kind=default_precision) :: zztop, dxx, dyy
logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
grids_mod::x_index
integer, parameter, public x_index
Definition: grids.F90:14
optionsdatabase_mod::options_get_integer
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
Definition: optionsdatabase.F90:217
grids_mod::y_index
integer, parameter, public y_index
Definition: grids.F90:14
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
grids_mod::z_index
integer, parameter, public z_index
Grid index parameters.
Definition: grids.F90:14
optionsdatabase_mod::options_get_logical
logical function, public options_get_logical(options_database, key, index)
Retrieves a logical value from the database that matches the provided key.
Definition: optionsdatabase.F90:154
optionsdatabase_mod::options_get_real
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
Definition: optionsdatabase.F90:91
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17