MONC
Data Types | Functions/Subroutines | Variables
casim_monc_dgs_space Module Reference

Data Types

type  casim_monc_dglist
 

Functions/Subroutines

subroutine allocate_casim_monc_dgs_space (current_state, casdiags)
 
subroutine populate_casim_monc_dg (current_state, casdiags)
 

Variables

type(casim_monc_dglistcasim_monc_dgs
 

Function/Subroutine Documentation

◆ allocate_casim_monc_dgs_space()

subroutine casim_monc_dgs_space::allocate_casim_monc_dgs_space ( type(model_state_type), intent(inout), target  current_state,
type(diaglist), intent(in), target  casdiags 
)

Definition at line 84 of file casim_monc_dgs_space.F90.

85 
86  type(model_state_type), target, intent(inout) :: current_state
87  type(diaglist), target, intent(in) :: casdiags
88 
89  integer :: z_size, y_size_local, x_size_local
90 
91  z_size = current_state%local_grid%size(z_index)
92  y_size_local = current_state%local_grid%size(y_index)
93  x_size_local = current_state%local_grid%size(x_index)
94 
95  ! precipitation arrays
96  if ( casdiags % l_precip ) then
97  allocate ( casim_monc_dgs % precip(y_size_local, x_size_local) )
98  casim_monc_dgs % precip(:,:) = 0.0_default_precision
99  endif
100 
101  if ( casdiags % l_surface_rain ) then
102  allocate ( casim_monc_dgs % SurfaceRainR(y_size_local, x_size_local) )
103  casim_monc_dgs % SurfaceRainR(:, :) = 0.0_default_precision
104  endif
105 
106  ! Process rate diagnostics
107  if ( casdiags % l_psedl ) then
108  allocate ( casim_monc_dgs % psedl(z_size, y_size_local, x_size_local) )
109  casim_monc_dgs % psedl(:, :, :) = 0.0_default_precision
110  endif
111 
112  if ( casdiags % l_pcond ) then
113  allocate ( casim_monc_dgs % pcond(z_size, y_size_local, x_size_local) )
114  casim_monc_dgs % pcond(:, :, :) = 0.0_default_precision
115  endif
116 
117  if ( casdiags % l_praut ) then
118  allocate ( casim_monc_dgs % praut(z_size, y_size_local, x_size_local) )
119  casim_monc_dgs % praut(:, :, :) = 0.0_default_precision
120  endif
121 
122  if ( casdiags % l_pracw ) then
123  allocate ( casim_monc_dgs % pracw(z_size, y_size_local, x_size_local) )
124  casim_monc_dgs % pracw(:, :, :) = 0.0_default_precision
125  endif
126 
127  if ( casdiags % l_prevp ) then
128  allocate ( casim_monc_dgs % prevp(z_size, y_size_local, x_size_local) )
129  casim_monc_dgs % prevp(:, :, :) = 0.0_default_precision
130  endif
131 
132  if ( casdiags % l_psedr ) then
133  allocate ( casim_monc_dgs % psedr(z_size, y_size_local, x_size_local) )
134  casim_monc_dgs % psedr(:, :, :) = 0.0_default_precision
135  endif
136 
137  if ( casdiags % l_dth ) then
138  ! potential temperature and mass tendencies
139  allocate ( casim_monc_dgs % dth_total(z_size, y_size_local, x_size_local) )
140  casim_monc_dgs % dth_total(:, :, :) = 0.0_default_precision
141 
142  allocate ( casim_monc_dgs % dth_cond_evap(z_size, y_size_local, x_size_local) )
143  casim_monc_dgs % dth_cond_evap(:, :, :) = 0.0_default_precision
144  endif
145 
146  if ( casdiags % l_dqv ) then
147  allocate ( casim_monc_dgs % dqv_total(z_size, y_size_local, x_size_local) )
148  casim_monc_dgs % dqv_total(:, :, :) = 0.0_default_precision
149 
150  allocate ( casim_monc_dgs % dqv_cond_evap(z_size, y_size_local, x_size_local) )
151  casim_monc_dgs % dqv_cond_evap(:, :, :) = 0.0_default_precision
152  endif
153 
154  if ( casdiags % l_dqc ) then
155  allocate ( casim_monc_dgs % dqc(z_size, y_size_local, x_size_local) )
156  casim_monc_dgs % dqc(:, :, :) = 0.0_default_precision
157  endif
158 
159  if ( casdiags % l_dqr ) then
160  allocate ( casim_monc_dgs % dqr(z_size, y_size_local, x_size_local) )
161  casim_monc_dgs % dqr(:, :, :) = 0.0_default_precision
162  endif
163 
164  if (.not. l_warm) then
165 
166  ! ice/cold process rates
167  if ( casdiags % l_surface_snow ) then
168  allocate ( casim_monc_dgs % SurfaceSnowR(y_size_local, x_size_local) )
169  casim_monc_dgs % SurfaceSnowR(:,:) = 0.0_default_precision
170  endif
171 
172  if ( casdiags % l_surface_graup ) then
173  allocate ( casim_monc_dgs % SurfaceGraupR(y_size_local, x_size_local) )
174  casim_monc_dgs % SurfaceGraupR(:,:) = 0.0_default_precision
175  endif
176 
177  if ( casdiags % l_phomc ) then
178  allocate ( casim_monc_dgs % phomc(z_size, y_size_local, x_size_local) )
179  casim_monc_dgs % phomc(:, :, :) = 0.0_default_precision
180  endif
181 
182  if ( casdiags % l_pinuc ) then
183  allocate ( casim_monc_dgs % pinuc(z_size, y_size_local, x_size_local) )
184  casim_monc_dgs % pinuc(:, :, :) = 0.0_default_precision
185  endif
186 
187  if ( casdiags % l_pidep ) then
188  allocate ( casim_monc_dgs % pidep(z_size, y_size_local, x_size_local) )
189  casim_monc_dgs % pidep(:, :, :) = 0.0_default_precision
190  endif
191 
192  if ( casdiags % l_piacw ) then
193  allocate ( casim_monc_dgs % piacw(z_size, y_size_local, x_size_local) )
194  casim_monc_dgs % piacw(:, :, :) = 0.0_default_precision
195  endif
196 
197  if ( casdiags % l_pisub ) then
198  allocate ( casim_monc_dgs % pisub(z_size, y_size_local, x_size_local) )
199  casim_monc_dgs % pisub(:, :, :) = 0.0_default_precision
200  endif
201 
202  if ( casdiags % l_pimlt ) then
203  allocate ( casim_monc_dgs % pimlt(z_size, y_size_local, x_size_local) )
204  casim_monc_dgs % pimlt(:, :, :) = 0.0_default_precision
205  endif
206 
207  if ( casdiags % l_psedi ) then
208  allocate ( casim_monc_dgs % psedi(z_size, y_size_local, x_size_local) )
209  casim_monc_dgs % psedi(:, :, :) = 0.0_default_precision
210  endif
211 
212  if ( casdiags % l_psdep ) then
213  allocate ( casim_monc_dgs % psdep(z_size, y_size_local, x_size_local) )
214  casim_monc_dgs % psdep(:, :, :) = 0.0_default_precision
215  endif
216 
217  if ( casdiags % l_psacw ) then
218  allocate ( casim_monc_dgs % psacw(z_size, y_size_local, x_size_local) )
219  casim_monc_dgs % psacw(:, :, :) = 0.0_default_precision
220  endif
221 
222  if ( casdiags % l_psacr ) then
223  allocate ( casim_monc_dgs % psacr(z_size, y_size_local, x_size_local) )
224  casim_monc_dgs % psacr(:, :, :) = 0.0_default_precision
225  endif
226 
227  if ( casdiags % l_pssub ) then
228  allocate ( casim_monc_dgs % pssub(z_size, y_size_local, x_size_local) )
229  casim_monc_dgs % pssub(:, :, :) = 0.0_default_precision
230  endif
231 
232  if ( casdiags % l_psmlt ) then
233  allocate ( casim_monc_dgs % psmlt(z_size, y_size_local, x_size_local) )
234  casim_monc_dgs % psmlt(:, :, :) = 0.0_default_precision
235  endif
236 
237  if ( casdiags % l_psaut ) then
238  allocate ( casim_monc_dgs % psaut(z_size, y_size_local, x_size_local) )
239  casim_monc_dgs % psaut(:, :, :) = 0.0_default_precision
240  endif
241 
242  if ( casdiags % l_psaci ) then
243  allocate ( casim_monc_dgs % psaci(z_size, y_size_local, x_size_local) )
244  casim_monc_dgs % psaci(:, :, :) = 0.0_default_precision
245  endif
246 
247  if ( casdiags % l_pgacw ) then
248  allocate ( casim_monc_dgs % pgacw(z_size, y_size_local, x_size_local) )
249  casim_monc_dgs % pgacw(:, :, :) = 0.0_default_precision
250  endif
251 
252  if ( casdiags % l_pgacs ) then
253  allocate ( casim_monc_dgs % pgacs(z_size, y_size_local, x_size_local) )
254  casim_monc_dgs % pgacs(:, :, :) = 0.0_default_precision
255  endif
256 
257  if ( casdiags % l_pgmlt ) then
258  allocate ( casim_monc_dgs % pgmlt(z_size, y_size_local, x_size_local) )
259  casim_monc_dgs % pgmlt(:, :, :) = 0.0_default_precision
260  endif
261 
262  if ( casdiags % l_pgsub ) then
263  allocate ( casim_monc_dgs % pgsub(z_size, y_size_local, x_size_local) )
264  casim_monc_dgs % pgsub(:, :, :) = 0.0_default_precision
265  endif
266 
267  if ( casdiags % l_pseds ) then
268  allocate ( casim_monc_dgs % pseds(z_size, y_size_local, x_size_local) )
269  casim_monc_dgs % pseds(:, :, :) = 0.0_default_precision
270  endif
271 
272  if ( casdiags % l_psedg ) then
273  allocate ( casim_monc_dgs % psedg(z_size, y_size_local, x_size_local) )
274  casim_monc_dgs % psedg(:, :, :) = 0.0_default_precision
275  endif
276 
277  ! ice, snow, graupel tendencies
278  if ( casdiags % l_dqi ) then
279  allocate ( casim_monc_dgs % dqi(z_size, y_size_local, x_size_local) )
280  casim_monc_dgs % dqi(:, :, :) = 0.0_default_precision
281  endif
282 
283  if ( casdiags % l_dqs ) then
284  allocate ( casim_monc_dgs % dqs(z_size, y_size_local, x_size_local) )
285  casim_monc_dgs % dqs(:, :, :) = 0.0_default_precision
286  endif
287 
288  if ( casdiags % l_dqg ) then
289  allocate ( casim_monc_dgs % dqg(z_size, y_size_local, x_size_local) )
290  casim_monc_dgs % dqg(:, :, :) = 0.0_default_precision
291  endif
292 
293  endif
294 

◆ populate_casim_monc_dg()

subroutine casim_monc_dgs_space::populate_casim_monc_dg ( type(model_state_type), intent(in), target  current_state,
type(diaglist), intent(in), target  casdiags 
)

Definition at line 297 of file casim_monc_dgs_space.F90.

298 
299  type(model_state_type), target, intent(in) :: current_state
300  type(diaglist), target, intent(in) :: casdiags
301 
302  INTEGER :: icol, jcol, target_x_index, target_y_index, k
303 
304  icol=current_state%column_local_x
305  jcol=current_state%column_local_y
306  target_y_index=jcol-current_state%local_grid%halo_size(y_index)
307  target_x_index=icol-current_state%local_grid%halo_size(x_index)
308 
309  if ( casdiags % l_surface_rain .and. casdiags % l_precip) &
310  casim_monc_dgs % precip(target_y_index,target_x_index) = &
311  casdiags % SurfaceRainR(1,1)
312 
313  if ( casdiags % l_surface_rain ) &
314  casim_monc_dgs % SurfaceRainR(target_y_index,target_x_index) = &
315  casdiags % SurfaceRainR(1,1)
316 
317  if ( casdiags % l_pcond ) &
318  casim_monc_dgs % pcond(:,target_y_index,target_x_index) = &
319  casdiags % pcond(1,1,:)
320  if ( casdiags % l_psedl ) &
321  casim_monc_dgs % psedl(:,target_y_index,target_x_index) = &
322  casdiags % psedl(1,1,:)
323  if ( casdiags % l_praut ) &
324  casim_monc_dgs % praut(:,target_y_index,target_x_index) = &
325  casdiags % praut(1,1,:)
326  if ( casdiags % l_pracw ) &
327  casim_monc_dgs % pracw(:,target_y_index,target_x_index) = &
328  casdiags % pracw(1,1,:)
329  if ( casdiags % l_prevp ) &
330  casim_monc_dgs % prevp(:,target_y_index,target_x_index) = &
331  casdiags % prevp(1,1,:)
332  if ( casdiags % l_psedr ) &
333  casim_monc_dgs % psedr(:,target_y_index,target_x_index) = &
334  casdiags % psedr(1,1,:)
335 
336  ! potential temperature and mass tendencies
337  if ( casdiags % l_dth ) then
338  casim_monc_dgs % dth_total(:,target_y_index,target_x_index) = &
339  casdiags % dth_total(1,1,:)
340  casim_monc_dgs % dth_cond_evap(:,target_y_index,target_x_index) = &
341  casdiags % dth_cond_evap(1,1,:)
342  endif
343  if ( casdiags % l_dqv ) then
344  casim_monc_dgs % dqv_total(:,target_y_index,target_x_index) = &
345  casdiags % dqv_total(1,1,:)
346  casim_monc_dgs % dqv_cond_evap(:,target_y_index,target_x_index) = &
347  casdiags % dqv_cond_evap(1,1,:)
348  endif
349  if ( casdiags % l_dqc ) &
350  casim_monc_dgs % dqc(:,target_y_index,target_x_index) = &
351  casdiags % dqc(1,1,:)
352  if ( casdiags % l_dqr ) &
353  casim_monc_dgs % dqr(:,target_y_index,target_x_index) = &
354  casdiags % dqr(1,1,:)
355 
356  if (.not. l_warm) then
357  if ( casdiags % l_precip ) &
358  casim_monc_dgs % precip(target_y_index,target_x_index) = &
359  casdiags % SurfaceRainR(1,1) + casdiags % SurfaceSnowR(1,1)
360  if ( casdiags % l_surface_snow ) &
361  casim_monc_dgs % SurfaceSnowR(target_y_index,target_x_index) = &
362  casdiags % SurfaceSnowR(1,1)
363  if ( casdiags % l_surface_graup ) &
364  casim_monc_dgs % SurfaceGraupR(target_y_index,target_x_index) = &
365  casdiags % SurfaceGraupR(1,1)
366  if ( casdiags % l_phomc ) &
367  casim_monc_dgs % phomc(:,target_y_index,target_x_index) = &
368  casdiags % psedr(1,1,:)
369  if ( casdiags % l_pinuc ) &
370  casim_monc_dgs % pinuc(:,target_y_index,target_x_index) = &
371  casdiags % pinuc(1,1,:)
372  if ( casdiags % l_pidep ) &
373  casim_monc_dgs % pidep(:,target_y_index,target_x_index) = &
374  casdiags % pidep(1,1,:)
375  if ( casdiags % l_piacw ) &
376  casim_monc_dgs % piacw(:,target_y_index,target_x_index) = &
377  casdiags % piacw(1,1,:)
378  if ( casdiags % l_pisub ) &
379  casim_monc_dgs % pisub(:,target_y_index,target_x_index) = &
380  casdiags % pisub(1,1,:)
381  if ( casdiags % l_pimlt ) &
382  casim_monc_dgs % pimlt(:,target_y_index,target_x_index) = &
383  casdiags % pimlt(1,1,:)
384  if ( casdiags % l_psedi ) &
385  casim_monc_dgs % psedi(:,target_y_index,target_x_index) = &
386  casdiags % psedi(1,1,:)
387  if ( casdiags % l_psmlt ) &
388  casim_monc_dgs % psmlt(:,target_y_index,target_x_index) = &
389  casdiags % psmlt(1,1,:)
390  if ( casdiags % l_psaut ) &
391  casim_monc_dgs % psaut(:,target_y_index,target_x_index) = &
392  casdiags % psaut(1,1,:)
393  if ( casdiags % l_psaci ) &
394  casim_monc_dgs % psaci(:,target_y_index,target_x_index) = &
395  casdiags % psaci(1,1,:)
396  if ( casdiags % l_psacw ) &
397  casim_monc_dgs % psacw(:,target_y_index,target_x_index) = &
398  casdiags % psacw(1,1,:)
399  if ( casdiags % l_psacr ) &
400  casim_monc_dgs % psacr(:,target_y_index,target_x_index) = &
401  casdiags % psacr(1,1,:)
402  if ( casdiags % l_pssub ) &
403  casim_monc_dgs % pssub(:,target_y_index,target_x_index) = &
404  casdiags % pssub(1,1,:)
405  if ( casdiags % l_psdep ) &
406  casim_monc_dgs % psdep(:,target_y_index,target_x_index) = &
407  casdiags % psdep(1,1,:)
408  if ( casdiags % l_pseds ) &
409  casim_monc_dgs % pseds(:,target_y_index,target_x_index) = &
410  casdiags % pseds(1,1,:)
411  if ( casdiags % l_pgacw ) &
412  casim_monc_dgs % pgacw(:,target_y_index,target_x_index) = &
413  casdiags % pgacw(1,1,:)
414  if ( casdiags % l_pgacs ) &
415  casim_monc_dgs % pgacs(:,target_y_index,target_x_index) = &
416  casdiags % pgacs(1,1,:)
417  if ( casdiags % l_pgmlt ) &
418  casim_monc_dgs % pgmlt(:,target_y_index,target_x_index) = &
419  casdiags % pgmlt(1,1,:)
420  if ( casdiags % l_pgsub ) &
421  casim_monc_dgs % pgsub(:,target_y_index,target_x_index) = &
422  casdiags % pgsub(1,1,:)
423  if ( casdiags % l_psedg ) &
424  casim_monc_dgs % psedg(:,target_y_index,target_x_index) = &
425  casdiags % psedg(1,1,:)
426  if ( casdiags % l_dqi ) &
427  casim_monc_dgs % dqi(:,target_y_index,target_x_index) = &
428  casdiags % dqi(1,1,:)
429  if ( casdiags % l_dqs ) &
430  casim_monc_dgs % dqs(:,target_y_index,target_x_index) = &
431  casdiags % dqs(1,1,:)
432  if ( casdiags % l_dqg ) &
433  casim_monc_dgs % dqg(:,target_y_index,target_x_index) = &
434  casdiags % dqg(1,1,:)
435 
436  endif
437 

Variable Documentation

◆ casim_monc_dgs

type (casim_monc_dglist) casim_monc_dgs_space::casim_monc_dgs

Definition at line 80 of file casim_monc_dgs_space.F90.

80  type (casim_monc_dglist) :: casim_monc_dgs