ROMS
Loading...
Searching...
No Matches
mod_coupling.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#ifdef SOLVE3D
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! DU_avg1 Time averaged U-flux for 2D equations (m3/s). !
13! DU_avg2 Time averaged U-flux for 3D equations coupling (m3/s). !
14! DV_avg1 Time averaged V-flux for 2D equations (m3/s). !
15! DV_avg2 Time averaged V-flux for 3D equations coupling (m3/s). !
16! Zt_avg1 Free-surface averaged over all short time-steps (m). !
17! rhoA Normalized vertical averaged density. !
18! rhoS Normalized vertical averaged density perturbation. !
19! rufrc Right-hand-side forcing term for 2D U-momentum (m4/s2) !
20! rvfrc Right-hand-side forcing term for 2D V-momentum (m4/s2) !
21! !
22!=======================================================================
23!
24 USE mod_kinds
25!
26 implicit none
27!
28 PUBLIC :: allocate_coupling
29 PUBLIC :: deallocate_coupling
30 PUBLIC :: initialize_coupling
31!
32!-----------------------------------------------------------------------
33! Define T_COUPLING structure.
34!-----------------------------------------------------------------------
35!
37!
38! Nonlinear model state.
39!
40 real(r8), pointer :: du_avg1(:,:)
41 real(r8), pointer :: du_avg2(:,:)
42 real(r8), pointer :: dv_avg1(:,:)
43 real(r8), pointer :: dv_avg2(:,:)
44 real(r8), pointer :: zt_avg1(:,:)
45 real(r8), pointer :: rufrc(:,:)
46 real(r8), pointer :: rvfrc(:,:)
47# ifdef VAR_RHO_2D
48 real(r8), pointer :: rhoa(:,:)
49 real(r8), pointer :: rhos(:,:)
50# endif
51
52# if defined TANGENT || defined TL_IOMS
53!
54! Tangent linear model state.
55!
56 real(r8), pointer :: tl_du_avg1(:,:)
57 real(r8), pointer :: tl_du_avg2(:,:)
58 real(r8), pointer :: tl_dv_avg1(:,:)
59 real(r8), pointer :: tl_dv_avg2(:,:)
60 real(r8), pointer :: tl_zt_avg1(:,:)
61 real(r8), pointer :: tl_rufrc(:,:)
62 real(r8), pointer :: tl_rvfrc(:,:)
63# ifdef VAR_RHO_2D_NOT_YET
64 real(r8), pointer :: tl_rhoa(:,:)
65 real(r8), pointer :: tl_rhos(:,:)
66# endif
67# endif
68
69# ifdef ADJOINT
70!
71! Adjoint model state.
72!
73 real(r8), pointer :: ad_du_avg1(:,:)
74 real(r8), pointer :: ad_du_avg2(:,:)
75 real(r8), pointer :: ad_dv_avg1(:,:)
76 real(r8), pointer :: ad_dv_avg2(:,:)
77 real(r8), pointer :: ad_zt_avg1(:,:)
78 real(r8), pointer :: ad_rufrc(:,:)
79 real(r8), pointer :: ad_rvfrc(:,:)
80# ifdef VAR_RHO_2D_NOT_YET
81 real(r8), pointer :: ad_rhoa(:,:)
82 real(r8), pointer :: ad_rhos(:,:)
83# endif
84# endif
85
86# if defined FORWARD_READ && \
87 (defined tangent || defined tl_ioms || defined adjoint)
88!
89! Latest two records of the nonlinear trajectory used to interpolate
90! the background state in the tangent linear and adjoint models.
91!
92 real(r8), pointer :: du_avg1g(:,:,:)
93 real(r8), pointer :: du_avg2g(:,:,:)
94 real(r8), pointer :: dv_avg1g(:,:,:)
95 real(r8), pointer :: dv_avg2g(:,:,:)
96 real(r8), pointer :: rufrcg(:,:,:)
97 real(r8), pointer :: rvfrcg(:,:,:)
98# endif
99
100 END TYPE t_coupling
101!
102 TYPE (t_coupling), allocatable :: coupling(:)
103!
104 CONTAINS
105!
106 SUBROUTINE allocate_coupling (ng, LBi, UBi, LBj, UBj)
107!
108!=======================================================================
109! !
110! This routine allocates all variables in the module for all nested !
111! grids. !
112! !
113!=======================================================================
114!
115 USE mod_param
116!
117! Imported variable declarations.
118!
119 integer, intent(in) :: ng, lbi, ubi, lbj, ubj
120!
121! Local variable declarations.
122!
123 real(r8) :: size2d
124!
125!-----------------------------------------------------------------------
126! Initialize module variables.
127!-----------------------------------------------------------------------
128!
129 IF (ng.eq.1) allocate ( coupling(ngrids) )
130!
131! Set horizontal array size.
132!
133 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
134!
135! Nonlinear model state.
136!
137 allocate ( coupling(ng) % DU_avg1(lbi:ubi,lbj:ubj) )
138 dmem(ng)=dmem(ng)+size2d
139
140 allocate ( coupling(ng) % DU_avg2(lbi:ubi,lbj:ubj) )
141 dmem(ng)=dmem(ng)+size2d
142
143 allocate ( coupling(ng) % DV_avg1(lbi:ubi,lbj:ubj) )
144 dmem(ng)=dmem(ng)+size2d
145
146 allocate ( coupling(ng) % DV_avg2(lbi:ubi,lbj:ubj) )
147 dmem(ng)=dmem(ng)+size2d
148
149 allocate ( coupling(ng) % Zt_avg1(lbi:ubi,lbj:ubj) )
150 dmem(ng)=dmem(ng)+size2d
151
152 allocate ( coupling(ng) % rufrc(lbi:ubi,lbj:ubj) )
153 dmem(ng)=dmem(ng)+size2d
154
155 allocate ( coupling(ng) % rvfrc(lbi:ubi,lbj:ubj) )
156 dmem(ng)=dmem(ng)+size2d
157
158# ifdef VAR_RHO_2D
159 allocate ( coupling(ng) % rhoA(lbi:ubi,lbj:ubj) )
160 dmem(ng)=dmem(ng)+size2d
161
162 allocate ( coupling(ng) % rhoS(lbi:ubi,lbj:ubj) )
163 dmem(ng)=dmem(ng)+size2d
164# endif
165
166# if defined TANGENT || defined TL_IOMS
167!
168! Tangent linear model state.
169!
170 allocate ( coupling(ng) % tl_DU_avg1(lbi:ubi,lbj:ubj) )
171 dmem(ng)=dmem(ng)+size2d
172
173 allocate ( coupling(ng) % tl_DU_avg2(lbi:ubi,lbj:ubj) )
174 dmem(ng)=dmem(ng)+size2d
175
176 allocate ( coupling(ng) % tl_DV_avg1(lbi:ubi,lbj:ubj) )
177 dmem(ng)=dmem(ng)+size2d
178
179 allocate ( coupling(ng) % tl_DV_avg2(lbi:ubi,lbj:ubj) )
180 dmem(ng)=dmem(ng)+size2d
181
182 allocate ( coupling(ng) % tl_Zt_avg1(lbi:ubi,lbj:ubj) )
183 dmem(ng)=dmem(ng)+size2d
184
185 allocate ( coupling(ng) % tl_rufrc(lbi:ubi,lbj:ubj) )
186 dmem(ng)=dmem(ng)+size2d
187
188 allocate ( coupling(ng) % tl_rvfrc(lbi:ubi,lbj:ubj) )
189 dmem(ng)=dmem(ng)+size2d
190
191# ifdef VAR_RHO_2D_NOT_YET
192 allocate ( coupling(ng) % tl_rhoA(lbi:ubi,lbj:ubj) )
193 dmem(ng)=dmem(ng)+size2d
194
195 allocate ( coupling(ng) % tl_rhoS(lbi:ubi,lbj:ubj) )
196 dmem(ng)=dmem(ng)+size2d
197# endif
198# endif
199
200# ifdef ADJOINT
201!
202! Adjoint model state.
203!
204 allocate ( coupling(ng) % ad_DU_avg1(lbi:ubi,lbj:ubj) )
205 dmem(ng)=dmem(ng)+size2d
206
207 allocate ( coupling(ng) % ad_DU_avg2(lbi:ubi,lbj:ubj) )
208 dmem(ng)=dmem(ng)+size2d
209
210 allocate ( coupling(ng) % ad_DV_avg1(lbi:ubi,lbj:ubj) )
211 dmem(ng)=dmem(ng)+size2d
212
213 allocate ( coupling(ng) % ad_DV_avg2(lbi:ubi,lbj:ubj) )
214 dmem(ng)=dmem(ng)+size2d
215
216 allocate ( coupling(ng) % ad_Zt_avg1(lbi:ubi,lbj:ubj) )
217 dmem(ng)=dmem(ng)+size2d
218
219 allocate ( coupling(ng) % ad_rufrc(lbi:ubi,lbj:ubj) )
220 dmem(ng)=dmem(ng)+size2d
221
222 allocate ( coupling(ng) % ad_rvfrc(lbi:ubi,lbj:ubj) )
223 dmem(ng)=dmem(ng)+size2d
224
225# ifdef VAR_RHO_2D_NOT_YET
226 allocate ( coupling(ng) % ad_rhoA(lbi:ubi,lbj:ubj) )
227 dmem(ng)=dmem(ng)+size2d
228
229 allocate ( coupling(ng) % ad_rhoS(lbi:ubi,lbj:ubj) )
230 dmem(ng)=dmem(ng)+size2d
231# endif
232# endif
233
234# if defined FORWARD_READ && \
235 (defined tangent || defined tl_ioms || defined adjoint)
236!
237! Latest two records of the nonlinear trajectory used to interpolate
238! the background state in the tangent linear and adjoint models.
239!
240 allocate ( coupling(ng) % DU_avg1G(lbi:ubi,lbj:ubj,2) )
241 dmem(ng)=dmem(ng)+2.0_r8*size2d
242
243 allocate ( coupling(ng) % DU_avg2G(lbi:ubi,lbj:ubj,2) )
244 dmem(ng)=dmem(ng)+2.0_r8*size2d
245
246 allocate ( coupling(ng) % DV_avg1G(lbi:ubi,lbj:ubj,2) )
247 dmem(ng)=dmem(ng)+2.0_r8*size2d
248
249 allocate ( coupling(ng) % DV_avg2G(lbi:ubi,lbj:ubj,2) )
250 dmem(ng)=dmem(ng)+2.0_r8*size2d
251
252 allocate ( coupling(ng) % rufrcG(lbi:ubi,lbj:ubj,2) )
253 dmem(ng)=dmem(ng)+2.0_r8*size2d
254
255 allocate ( coupling(ng) % rvfrcG(lbi:ubi,lbj:ubj,2) )
256 dmem(ng)=dmem(ng)+2.0_r8*size2d
257# endif
258!
259 RETURN
260 END SUBROUTINE allocate_coupling
261!
262 SUBROUTINE deallocate_coupling (ng)
263!
264!=======================================================================
265! !
266! This routine deallocates all variables in the module for all nested !
267! grids. !
268! !
269!=======================================================================
270!
271 USE mod_param, ONLY : ngrids
272# ifdef SUBOBJECT_DEALLOCATION
273 USE destroy_mod, ONLY : destroy
274# endif
275!
276! Imported variable declarations.
277!
278 integer, intent(in) :: ng
279!
280! Local variable declarations.
281!
282 character (len=*), parameter :: myfile = &
283 & __FILE__//", deallocate_coupling"
284
285# ifdef SUBOBJECT_DEALLOCATION
286!
287!-----------------------------------------------------------------------
288! Deallocate module variables.
289!-----------------------------------------------------------------------
290!
291! Nonlinear model state.
292!
293 IF (.not.destroy(ng, coupling(ng)%DU_avg1, myfile, &
294 & __line__, 'COUPLING(ng)%DU_avg1')) RETURN
295
296 IF (.not.destroy(ng, coupling(ng)%DU_avg2, myfile, &
297 & __line__, 'COUPLING(ng)%DU_avg2')) RETURN
298
299 IF (.not.destroy(ng, coupling(ng)%DV_avg1, myfile, &
300 & __line__, 'COUPLING(ng)%DV_avg1')) RETURN
301
302 IF (.not.destroy(ng, coupling(ng)%DV_avg2, myfile, &
303 & __line__, 'COUPLING(ng)%DV_avg2')) RETURN
304
305 IF (.not.destroy(ng, coupling(ng)%Zt_avg1, myfile, &
306 & __line__, 'COUPLING(ng)%Zt_avg1')) RETURN
307
308 IF (.not.destroy(ng, coupling(ng)%rufrc, myfile, &
309 & __line__, 'COUPLING(ng)%rufrc')) RETURN
310
311 IF (.not.destroy(ng, coupling(ng)%rvfrc, myfile, &
312 & __line__, 'COUPLING(ng)%rvfrc')) RETURN
313
314# ifdef VAR_RHO_2D
315 IF (.not.destroy(ng, coupling(ng)%rhoA, myfile, &
316 & __line__, 'COUPLING(ng)%rhoA')) RETURN
317
318 IF (.not.destroy(ng, coupling(ng)%rhoS, myfile, &
319 & __line__, 'COUPLING(ng)%rhoS')) RETURN
320# endif
321
322# if defined TANGENT || defined TL_IOMS
323!
324! Tangent linear model state.
325!
326 IF (.not.destroy(ng, coupling(ng)%tl_DU_avg1, myfile, &
327 & __line__, 'COUPLING(ng)%tl_DU_avg1')) RETURN
328
329 IF (.not.destroy(ng, coupling(ng)%tl_DU_avg2, myfile, &
330 & __line__, 'COUPLING(ng)%tl_DU_avg2')) RETURN
331
332 IF (.not.destroy(ng, coupling(ng)%tl_DV_avg1, myfile, &
333 & __line__, 'COUPLING(ng)%tl_DV_avg1')) RETURN
334
335 IF (.not.destroy(ng, coupling(ng)%tl_DV_avg2, myfile, &
336 & __line__, 'COUPLING(ng)%tl_DV_avg2')) RETURN
337
338 IF (.not.destroy(ng, coupling(ng)%tl_Zt_avg1, myfile, &
339 & __line__, 'COUPLING(ng)%tl_Zt_avg1')) RETURN
340
341 IF (.not.destroy(ng, coupling(ng)%tl_rufrc, myfile, &
342 & __line__, 'COUPLING(ng)%tl_rufrc')) RETURN
343
344 IF (.not.destroy(ng, coupling(ng)%tl_rvfrc, myfile, &
345 & __line__, 'COUPLING(ng)%tl_rvfrc')) RETURN
346
347# ifdef VAR_RHO_2D_NOT_YET
348 IF (.not.destroy(ng, coupling(ng)%tl_rhoA, myfile, &
349 & __line__, 'COUPLING(ng)%tl_rhoA')) RETURN
350
351 IF (.not.destroy(ng, coupling(ng)%tl_rhoS, myfile, &
352 & __line__, 'COUPLING(ng)%tl_rhoS')) RETURN
353# endif
354# endif
355
356# ifdef ADJOINT
357!
358! Adjoint model state.
359!
360 IF (.not.destroy(ng, coupling(ng)%ad_DU_avg1, myfile, &
361 & __line__, 'COUPLING(ng)%ad_DU_avg1')) RETURN
362
363 IF (.not.destroy(ng, coupling(ng)%ad_DU_avg2, myfile, &
364 & __line__, 'COUPLING(ng)%ad_DU_avg2')) RETURN
365
366 IF (.not.destroy(ng, coupling(ng)%ad_DV_avg1, myfile, &
367 & __line__, 'COUPLING(ng)%ad_DV_avg1')) RETURN
368
369 IF (.not.destroy(ng, coupling(ng)%ad_DV_avg2, myfile, &
370 & __line__, 'COUPLING(ng)%ad_DV_avg2')) RETURN
371
372 IF (.not.destroy(ng, coupling(ng)%ad_Zt_avg1, myfile, &
373 & __line__, 'COUPLING(ng)%ad_Zt_avg1')) RETURN
374
375 IF (.not.destroy(ng, coupling(ng)%ad_rufrc, myfile, &
376 & __line__, 'COUPLING(ng)%ad_rufrc')) RETURN
377
378 IF (.not.destroy(ng, coupling(ng)%ad_rvfrc, myfile, &
379 & __line__, 'COUPLING(ng)%ad_rvfrc')) RETURN
380
381# ifdef VAR_RHO_2D_NOT_YET
382 IF (.not.destroy(ng, coupling(ng)%ad_rhoA, myfile, &
383 & __line__, 'COUPLING(ng)%ad_rhoA')) RETURN
384
385 IF (.not.destroy(ng, coupling(ng)%ad_rhoS, myfile, &
386 & __line__, 'COUPLING(ng)%ad_rhoS')) RETURN
387# endif
388# endif
389
390# if defined FORWARD_READ && \
391 (defined tangent || defined tl_ioms || defined adjoint)
392!
393! Latest two records of the nonlinear trajectory used to interpolate
394! the background state in the tangent linear and adjoint models.
395!
396 IF (.not.destroy(ng, coupling(ng)%DU_avg1G, myfile, &
397 & __line__, 'COUPLING(ng)%DU_avg1G')) RETURN
398
399 IF (.not.destroy(ng, coupling(ng)%DU_avg2G, myfile, &
400 & __line__, 'COUPLING(ng)%DU_avg2G')) RETURN
401
402 IF (.not.destroy(ng, coupling(ng)%DV_avg1G, myfile, &
403 & __line__, 'COUPLING(ng)%DV_avg1G')) RETURN
404
405 IF (.not.destroy(ng, coupling(ng)%DV_avg2G, myfile, &
406 & __line__, 'COUPLING(ng)%DV_avg2G')) RETURN
407
408 IF (.not.destroy(ng, coupling(ng)%rufrcG, myfile, &
409 & __line__, 'COUPLING(ng)%rufrcG')) RETURN
410
411 IF (.not.destroy(ng, coupling(ng)%rvfrcG, myfile, &
412 & __line__, 'COUPLING(ng)%rvfrcG')) RETURN
413# endif
414# endif
415!
416!-----------------------------------------------------------------------
417! Deallocate derived-type COUPLING structure.
418!-----------------------------------------------------------------------
419!
420 IF (ng.eq.ngrids) THEN
421 IF (allocated(coupling)) deallocate ( coupling )
422 END IF
423!
424 RETURN
425 END SUBROUTINE deallocate_coupling
426!
427 SUBROUTINE initialize_coupling (ng, tile, model)
428!
429!=======================================================================
430! !
431! This routine initialize all variables in the module using first !
432! touch distribution policy. In shared-memory configuration, this !
433! operation actually performs propagation of the "shared arrays" !
434! across the cluster, unless another policy is specified to !
435! override the default. !
436! !
437!=======================================================================
438!
439 USE mod_param
440!
441! Imported variable declarations.
442!
443 integer, intent(in) :: ng, tile, model
444!
445! Local variable declarations.
446!
447 integer :: imin, imax, jmin, jmax
448 integer :: i, j
449
450 real(r8), parameter :: inival = 0.0_r8
451
452# include "set_bounds.h"
453!
454! Set array initialization range.
455!
456# ifdef DISTRIBUTE
457 imin=bounds(ng)%LBi(tile)
458 imax=bounds(ng)%UBi(tile)
459 jmin=bounds(ng)%LBj(tile)
460 jmax=bounds(ng)%UBj(tile)
461# else
462 IF (domain(ng)%Western_Edge(tile)) THEN
463 imin=bounds(ng)%LBi(tile)
464 ELSE
465 imin=istr
466 END IF
467 IF (domain(ng)%Eastern_Edge(tile)) THEN
468 imax=bounds(ng)%UBi(tile)
469 ELSE
470 imax=iend
471 END IF
472 IF (domain(ng)%Southern_Edge(tile)) THEN
473 jmin=bounds(ng)%LBj(tile)
474 ELSE
475 jmin=jstr
476 END IF
477 IF (domain(ng)%Northern_Edge(tile)) THEN
478 jmax=bounds(ng)%UBj(tile)
479 ELSE
480 jmax=jend
481 END IF
482# endif
483!
484!-----------------------------------------------------------------------
485! Initialize module variables.
486!-----------------------------------------------------------------------
487!
488! Nonlinear model state.
489!
490 IF ((model.eq.0).or.(model.eq.inlm)) THEN
491 DO j=jmin,jmax
492 DO i=imin,imax
493 coupling(ng) % DU_avg1(i,j) = inival
494 coupling(ng) % DU_avg2(i,j) = inival
495
496 coupling(ng) % DV_avg1(i,j) = inival
497 coupling(ng) % DV_avg2(i,j) = inival
498
499 coupling(ng) % Zt_avg1(i,j) = inival
500
501 coupling(ng) % rufrc(i,j) = inival
502 coupling(ng) % rvfrc(i,j) = inival
503
504# ifdef VAR_RHO_2D
505 coupling(ng) % rhoA(i,j) = inival
506 coupling(ng) % rhoS(i,j) = inival
507# endif
508 END DO
509 END DO
510 END IF
511
512# if defined TANGENT || defined TL_IOMS
513!
514! Tangent linear model state.
515!
516 IF ((model.eq.0).or.(model.eq.itlm).or.(model.eq.irpm)) THEN
517 DO j=jmin,jmax
518 DO i=imin,imax
519 coupling(ng) % tl_DU_avg1(i,j) = inival
520 coupling(ng) % tl_DU_avg2(i,j) = inival
521
522 coupling(ng) % tl_DV_avg1(i,j) = inival
523 coupling(ng) % tl_DV_avg2(i,j) = inival
524
525 coupling(ng) % tl_Zt_avg1(i,j) = inival
526
527 coupling(ng) % tl_rufrc(i,j) = inival
528 coupling(ng) % tl_rvfrc(i,j) = inival
529
530# ifdef VAR_RHO_2D_NOT_YET
531 coupling(ng) % tl_rhoA(i,j) = inival
532 coupling(ng) % tl_rhoS(i,j) = inival
533# endif
534 END DO
535 END DO
536 END IF
537# endif
538
539# ifdef ADJOINT
540!
541! Adjoint model state.
542!
543 IF ((model.eq.0).or.(model.eq.iadm)) THEN
544 DO j=jmin,jmax
545 DO i=imin,imax
546 coupling(ng) % ad_DU_avg1(i,j) = inival
547 coupling(ng) % ad_DU_avg2(i,j) = inival
548
549 coupling(ng) % ad_DV_avg1(i,j) = inival
550 coupling(ng) % ad_DV_avg2(i,j) = inival
551
552 coupling(ng) % ad_Zt_avg1(i,j) = inival
553
554 coupling(ng) % ad_rufrc(i,j) = inival
555 coupling(ng) % ad_rvfrc(i,j) = inival
556
557# ifdef VAR_RHO_2D_NOT_YET
558 coupling(ng) % ad_rhoA(i,j) = inival
559 coupling(ng) % ad_rhoS(i,j) = inival
560# endif
561 END DO
562 END DO
563 END IF
564# endif
565
566# if defined FORWARD_READ && \
567 (defined tangent || defined tl_ioms || defined adjoint)
568!
569! Latest two records of the nonlinear trajectory used to interpolate
570! the background state in the tangent linear and adjoint models.
571!
572 IF (model.eq.0) THEN
573 DO j=jmin,jmax
574 DO i=imin,imax
575 coupling(ng) % DU_avg1G(i,j,1) = inival
576 coupling(ng) % DU_avg1G(i,j,2) = inival
577 coupling(ng) % DU_avg2G(i,j,1) = inival
578 coupling(ng) % DU_avg2G(i,j,2) = inival
579
580 coupling(ng) % DV_avg1G(i,j,1) = inival
581 coupling(ng) % DV_avg1G(i,j,2) = inival
582 coupling(ng) % DV_avg2G(i,j,1) = inival
583 coupling(ng) % DV_avg2G(i,j,2) = inival
584
585 coupling(ng) % rufrcG(i,j,1) = inival
586 coupling(ng) % rufrcG(i,j,2) = inival
587 coupling(ng) % rvfrcG(i,j,1) = inival
588 coupling(ng) % rvfrcG(i,j,2) = inival
589 END DO
590 END DO
591 END IF
592# endif
593!
594 RETURN
595 END SUBROUTINE initialize_coupling
596#endif
597 END MODULE mod_coupling
subroutine, public deallocate_coupling(ng)
type(t_coupling), dimension(:), allocatable coupling
subroutine, public initialize_coupling(ng, tile, model)
subroutine, public allocate_coupling(ng, lbi, ubi, lbj, ubj)
integer, parameter r8
Definition mod_kinds.F:28
integer, parameter inlm
Definition mod_param.F:662
integer, parameter irpm
Definition mod_param.F:664
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
integer, parameter iadm
Definition mod_param.F:665
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ngrids
Definition mod_param.F:113
integer, parameter itlm
Definition mod_param.F:663