ROMS
Loading...
Searching...
No Matches
mod_arrays.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 MODULE mod_arrays
3!
4!git $Id$
5!================================================== Hernan G. Arango ===
6! Copyright (c) 2002-2025 The ROMS Group !
7! Licensed under a MIT/X style license !
8! See License_ROMS.md !
9!=======================================================================
10! !
11! This module is used to allocate, initialize, and deallocate ROMS !
12! state arrays for each nested and/or multiple connected grids. !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_parallel
18 USE mod_scalars
19!
20#if defined AVERAGES || \
21 (defined ad_averages && defined adjoint) || \
22 (defined rp_averages && defined tl_ioms) || \
23 (defined tl_averages && defined tangent)
24 USE mod_average, ONLY : allocate_average, &
27#endif
28 USE mod_boundary, ONLY : allocate_boundary, &
31 USE mod_clima, ONLY : allocate_clima, &
34#ifdef SOLVE3D
35 USE mod_coupling, ONLY : allocate_coupling, &
38#endif
39#ifdef DIAGNOSTICS
40 USE mod_diags, ONLY : allocate_diags, &
43#endif
44#ifdef GRID_EXTRACT
45 USE mod_extract, ONLY : allocate_extract, &
46 & deallocate_extract, &
47 & initialize_extract
48#endif
49 USE mod_forces, ONLY : allocate_forces, &
52#if defined FOUR_DVAR || defined VERIFICATION
55#endif
56 USE mod_grid, ONLY : allocate_grid, &
59#ifdef ICE_MODEL
60 USE mod_ice, ONLY : allocate_ice, &
63#endif
65 USE mod_mixing, ONLY : allocate_mixing, &
68#ifdef NESTING
69 USE mod_nesting, ONLY : allocate_nesting, &
72#endif
73 USE mod_ocean, ONLY : allocate_ocean, &
76#if defined SEDIMENT || defined BBL_MODEL
77 USE mod_sedbed, ONLY : allocate_sedbed, &
80#endif
81 USE mod_sources, ONLY : allocate_sources, &
83#if defined SSH_TIDES || defined UV_TIDES
84 USE mod_tides, ONLY : allocate_tides, &
87#endif
88#ifdef BBL_MODEL
89 USE mod_bbl, ONLY : allocate_bbl, &
92#endif
93#if defined PIO_LIB && defined DISTRIBUTE
94 USE set_pio_mod, ONLY : set_iodecomp
95#endif
96!
97 implicit none
98!
99 PUBLIC :: roms_allocate_arrays
100 PUBLIC :: roms_deallocate_arrays
101 PUBLIC :: roms_initialize_arrays
102!
103#if defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
104 defined opt_observations || defined sensitivity_4dvar || \
105 defined so_semi
106 logical :: lallocateclima = .true.
107#else
108 logical :: lallocateclima = .false.
109#endif
110!
111 CONTAINS
112!
113 SUBROUTINE roms_allocate_arrays (allocate_vars)
114!
115!=======================================================================
116! !
117! This routine allocates ROMS state variables. !
118! !
119!=======================================================================
120!
121! Imported variable declarations
122!
123 logical, intent(in) :: allocate_vars
124!
125! Local variable declarations.
126!
127 integer :: ng, thread, tile
128 integer :: imins, imaxs, jmins, jmaxs
129 integer :: lbi, ubi, lbj, ubj, lbij, ubij
130!
131 character (len=*), parameter :: myfile = &
132 & __FILE__//", ROMS_allocate_arrays"
133
134#ifdef PROFILE
135!
136!-----------------------------------------------------------------------
137! Turn on allocation time wall clock.
138!-----------------------------------------------------------------------
139!
140 DO ng=1,ngrids
141 DO thread=thread_range
142 CALL wclock_on (ng, inlm, 1, __line__, myfile)
143 END DO
144!$OMP BARRIER
145 END DO
146#endif
147!
148!-----------------------------------------------------------------------
149! Allocate model type-derived structures and its variables.
150!-----------------------------------------------------------------------
151!
152 IF (allocate_vars) then
153#ifdef DISTRIBUTE
154 tile=myrank
155#else
156 tile=0
157#endif
158 DO ng=1,ngrids
159!$OMP MASTER
160 lbi=bounds(ng)%LBi(tile)
161 ubi=bounds(ng)%UBi(tile)
162 lbj=bounds(ng)%LBj(tile)
163 ubj=bounds(ng)%UBj(tile)
164 lbij=bounds(ng)%LBij
165 ubij=bounds(ng)%UBij
166#if defined AVERAGES || \
167 (defined ad_averages && defined adjoint) || \
168 (defined rp_averages && defined tl_ioms) || \
169 (defined tl_averages && defined tangent)
170 CALL allocate_average (ng, lbi, ubi, lbj, ubj)
171#endif
172 CALL allocate_boundary (ng)
173#ifdef BBL_MODEL
174 CALL allocate_bbl (ng, lbi, ubi, lbj, ubj)
175#endif
176 IF (lallocateclima.or.lclimatology(ng)) THEN
177 CALL allocate_clima (ng, lbi, ubi, lbj, ubj)
178 END IF
179#ifdef SOLVE3D
180 CALL allocate_coupling (ng, lbi, ubi, lbj, ubj)
181#endif
182#ifdef DIAGNOSTICS
183 CALL allocate_diags (ng, lbi, ubi, lbj, ubj)
184#endif
185#ifdef GRID_EXTRACT
186 CALL allocate_extract (ng, extractflag(ng), &
187 & xtr_bounds(ng)%LBi(tile), &
188 & xtr_bounds(ng)%UBi(tile), &
189 & xtr_bounds(ng)%LBj(tile), &
190 & xtr_bounds(ng)%UBj(tile))
191#endif
192 CALL allocate_forces (ng, lbi, ubi, lbj, ubj)
193 CALL allocate_grid (ng, extractflag(ng), &
194 & lbi, ubi, lbj, ubj, lbij, ubij)
195#ifdef ICE_MODEL
196 CALL allocate_ice (ng, lbi, ubi, lbj, ubj, .true.)
197#endif
198 CALL allocate_mixing (ng, lbi, ubi, lbj, ubj)
199 CALL allocate_ocean (ng, lbi, ubi, lbj, ubj)
200#if defined SEDIMENT || defined BBL_MODEL
201 CALL allocate_sedbed (ng, lbi, ubi, lbj, ubj)
202#endif
203#if defined SSH_TIDES || defined UV_TIDES
204 CALL allocate_tides (ng, lbi, ubi, lbj, ubj)
205#endif
206 IF (luvsrc(ng).or.lwsrc(ng).or.any(ltracersrc(:,ng))) THEN
207 CALL allocate_sources (ng)
208 END IF
209!$OMP END MASTER
210!$OMP BARRIER
211 END DO
212
213#ifdef NESTING
214!
215! Allocate and initialized contact points boundaty structure. It
216! needs to be delayed to the end because we need "LBC_apply" to
217! allocated in "mod_boundary" for all nested grid.
218!
219!$OMP MASTER
221!$OMP END MASTER
222!$OMP BARRIER
223#endif
224
225#if defined PIO_LIB && defined DISTRIBUTE
226!
227! Allocate and initialize PIO single and double precision data
228! decomposition for mapping between computational and I/O processes.
229! It needs to done after all parameters have been allocated and
230! initialize. In particular, we need tidal constituents and state
231! propagator parameters.
232!
233 CALL set_iodecomp
234#endif
235
236 lallocatedmemory=.true.
237
238 END IF
239
240#ifdef PROFILE
241!
242!-----------------------------------------------------------------------
243! Turn off allocation time wall clock.
244!-----------------------------------------------------------------------
245!
246 DO ng=1,ngrids
247 DO thread=thread_range
248 CALL wclock_off (ng, inlm, 1, __line__, myfile)
249 END DO
250!$OMP BARRIER
251 END DO
252#endif
253!
254 RETURN
255 END SUBROUTINE roms_allocate_arrays
256!
258!
259!=======================================================================
260! !
261! This routine deallocates ROMS state objects, arrays and vectors. !
262! !
263!=======================================================================
264!
265 USE mod_param, ONLY : ngrids
266!
267! Local variable declarations.
268!
269 integer :: ng
270!
271 character (len=*), parameter :: myfile = &
272 & __FILE__//", ROMS_deallocate_arrays"
273!
274!-----------------------------------------------------------------------
275! Deallocate all structures.
276!-----------------------------------------------------------------------
277!
278 DO ng=1,ngrids
279!$OMP MASTER
280#if defined AVERAGES || \
281 (defined ad_averages && defined adjoint) || \
282 (defined rp_averages && defined tl_ioms) || \
283 (defined tl_averages && defined tangent)
284 CALL deallocate_average (ng)
285#endif
286 CALL deallocate_boundary (ng)
287#ifdef BBL_MODEL
288 CALL deallocate_bbl (ng)
289#endif
290 IF (lallocateclima.or.lclimatology(ng)) THEN
291 CALL deallocate_clima (ng)
292 END IF
293#ifdef SOLVE3D
294 CALL deallocate_coupling (ng)
295#endif
296#ifdef DIAGNOSTICS
297 CALL deallocate_diags (ng)
298#endif
299#ifdef GRID_EXTRACT
300 CALL deallocate_extract (ng)
301#endif
302 CALL deallocate_forces (ng)
303 CALL deallocate_grid (ng)
304#ifdef ICE_MODEL
305 CALL deallocate_ice (ng)
306#endif
307 CALL deallocate_mixing (ng)
308 CALL deallocate_ocean (ng)
309#if defined SEDIMENT || defined BBL_MODEL
310 CALL deallocate_sedbed (ng)
311#endif
312#if defined SSH_TIDES || defined UV_TIDES
313 CALL deallocate_tides (ng)
314#endif
315 IF (luvsrc(ng).or.lwsrc(ng).or.any(ltracersrc(:,ng))) THEN
316 CALL deallocate_sources (ng)
317 END IF
318!$OMP END MASTER
319!$OMP BARRIER
320 END DO
321
322#ifdef NESTING
323!
324! Deallocate nesting variables and structures.
325!
326!$OMP MASTER
328!$OMP END MASTER
329!$OMP BARRIER
330#endif
331
332#if defined FOUR_DVAR || defined VERIFICATION
333!
334! Deallocate observations arrays and object.
335!
337#endif
338!
339! Deallocate I/O derived-type structures.
340!
342!
343! Deallocate main configuration dimensions and associated parameters.
344!
346!
347 RETURN
348 END SUBROUTINE roms_deallocate_arrays
349!
351!
352!=======================================================================
353! !
354! This routine initialize ROMS state variables. In shared-memory it !
355! important for the first-touch policy in memory. !
356! !
357!=======================================================================
358!
359! Local variable declarations.
360!
361 integer :: ng, thread, tile
362!
363 integer, parameter :: model = 0
364!
365 character (len=*), parameter :: myfile = &
366 & __FILE__//", ROMS_initialize_arrays"
367
368#ifdef PROFILE
369!
370!-----------------------------------------------------------------------
371! Turn on allocation time wall clock.
372!-----------------------------------------------------------------------
373!
374 DO ng=1,ngrids
375 DO thread=thread_range
376 CALL wclock_on (ng, inlm, 1, __line__, myfile)
377 END DO
378!$OMP BARRIER
379 END DO
380#endif
381!
382!-----------------------------------------------------------------------
383! Intialize variables within structures for each grid.
384!-----------------------------------------------------------------------
385!
386 DO ng=1,ngrids
387#ifdef NESTING
388 IF (ng.eq.1) THEN
390 END IF
391#endif
392 DO tile=first_tile(ng),last_tile(ng),+1
393#if defined AVERAGES || \
394 (defined ad_averages && defined adjoint) || \
395 (defined rp_averages && defined tl_ioms) || \
396 (defined tl_averages && defined tangent)
397 CALL initialize_average (ng, tile)
398#endif
399#ifdef BBL_MODEL
400 CALL initialize_bbl (ng, tile)
401#endif
402 CALL initialize_boundary (ng, tile, model)
403
404 IF (lallocateclima.or.lclimatology(ng)) THEN
405 CALL initialize_clima (ng, tile)
406 END IF
407#ifdef SOLVE3D
408 CALL initialize_coupling (ng, tile, model)
409#endif
410#ifdef DIAGNOSTICS
411 CALL initialize_diags (ng, tile)
412#endif
413#ifdef GRID_EXTRACT
414 CALL initialize_extract (ng, tile, model)
415#endif
416 CALL initialize_forces (ng, tile, model)
417 CALL initialize_grid (ng, tile, model)
418#ifdef ICE_MODEL
419 CALL initialize_ice (ng, tile, model)
420#endif
421 CALL initialize_mixing (ng, tile, model)
422 CALL initialize_ocean (ng, tile, model)
423#if defined SEDIMENT || defined BBL_MODEL
424 CALL initialize_sedbed (ng, tile, model)
425#endif
426#if defined SSH_TIDES || defined UV_TIDES
427 CALL initialize_tides (ng, tile)
428#endif
429 END DO
430!$OMP BARRIER
431 END DO
432
433#if defined FOUR_DVAR || defined VERIFICATION
434!
435! Finish allocating observation arrays. Then, initialize observation
436! arrays. Notice that some of the module variables are allocated in
437! the call to "allocate_fourdvar" in "read_phypar".
438!
440#endif
441
442#ifdef PROFILE
443!
444!-----------------------------------------------------------------------
445! Turn off allocation time wall clock.
446!-----------------------------------------------------------------------
447!
448 DO ng=1,ngrids
449 DO thread=thread_range
450 CALL wclock_off (ng, inlm, 1, __line__, myfile)
451 END DO
452!$OMP BARRIER
453 END DO
454#endif
455!
456 RETURN
457 END SUBROUTINE roms_initialize_arrays
458!
459 END MODULE mod_arrays
subroutine, public roms_initialize_arrays
Definition mod_arrays.F:351
subroutine, public roms_allocate_arrays(allocate_vars)
Definition mod_arrays.F:114
logical lallocateclima
Definition mod_arrays.F:106
subroutine, public roms_deallocate_arrays
Definition mod_arrays.F:258
subroutine, public deallocate_average(ng)
subroutine, public initialize_average(ng, tile)
subroutine, public allocate_average(ng, lbi, ubi, lbj, ubj)
subroutine, public allocate_bbl(ng, lbi, ubi, lbj, ubj)
Definition mod_bbl.F:67
subroutine, public initialize_bbl(ng, tile)
Definition mod_bbl.F:209
subroutine, public deallocate_bbl(ng)
Definition mod_bbl.F:134
subroutine, public initialize_boundary(ng, tile, model)
subroutine, public deallocate_boundary(ng)
subroutine, public allocate_boundary(ng)
subroutine, public allocate_clima(ng, lbi, ubi, lbj, ubj)
Definition mod_clima.F:158
subroutine, public deallocate_clima(ng)
Definition mod_clima.F:320
subroutine, public initialize_clima(ng, tile)
Definition mod_clima.F:492
subroutine, public deallocate_coupling(ng)
subroutine, public initialize_coupling(ng, tile, model)
subroutine, public allocate_coupling(ng, lbi, ubi, lbj, ubj)
subroutine, public initialize_diags(ng, tile)
Definition mod_diags.F:364
subroutine, public allocate_diags(ng, lbi, ubi, lbj, ubj)
Definition mod_diags.F:105
subroutine, public deallocate_diags(ng)
Definition mod_diags.F:231
subroutine, public initialize_forces(ng, tile, model)
subroutine, public deallocate_forces(ng)
subroutine, public allocate_forces(ng, lbi, ubi, lbj, ubj)
Definition mod_forces.F:559
subroutine, public deallocate_fourdvar
subroutine, public initialize_fourdvar
subroutine, public initialize_grid(ng, tile, model)
Definition mod_grid.F:1241
subroutine, public allocate_grid(ng, extract_flag, lbi, ubi, lbj, ubj, lbij, ubij)
Definition mod_grid.F:371
subroutine, public deallocate_grid(ng)
Definition mod_grid.F:838
subroutine, public deallocate_ice(ng)
Definition ice_mod.h:533
subroutine, public initialize_ice(ng, tile, model)
Definition ice_mod.h:647
subroutine, public allocate_ice(ng, lbi, ubi, lbj, ubj, ice_kernel)
Definition ice_mod.h:326
subroutine, public deallocate_iounits
subroutine, public initialize_mixing(ng, tile, model)
subroutine, public deallocate_mixing(ng)
Definition mod_mixing.F:832
subroutine, public allocate_mixing(ng, lbi, ubi, lbj, ubj)
Definition mod_mixing.F:404
subroutine, public deallocate_nesting
subroutine, public initialize_nesting
subroutine, public allocate_nesting
subroutine, public deallocate_ocean(ng)
Definition mod_ocean.F:941
subroutine, public initialize_ocean(ng, tile, model)
Definition mod_ocean.F:1526
subroutine, public allocate_ocean(ng, lbi, ubi, lbj, ubj)
Definition mod_ocean.F:356
integer, dimension(:), allocatable first_tile
integer, dimension(:), allocatable last_tile
integer, parameter inlm
Definition mod_param.F:662
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
subroutine, public deallocate_param
Definition mod_param.F:892
integer ngrids
Definition mod_param.F:113
logical, dimension(:), allocatable luvsrc
logical, dimension(:,:), allocatable ltracersrc
logical lallocatedmemory
logical, dimension(:), allocatable lwsrc
integer, dimension(:), allocatable extractflag
logical, dimension(:), allocatable lclimatology
subroutine, public deallocate_sedbed(ng)
Definition sedbed_mod.h:269
subroutine, public initialize_sedbed(ng, tile, model)
Definition sedbed_mod.h:445
subroutine, public allocate_sedbed(ng, lbi, ubi, lbj, ubj)
Definition sedbed_mod.h:162
subroutine, public allocate_sources(ng)
subroutine, public deallocate_sources(ng)
subroutine, public deallocate_tides(ng)
Definition mod_tides.F:345
subroutine, public initialize_tides(ng, tile)
Definition mod_tides.F:499
subroutine, public allocate_tides(ng, lbi, ubi, lbj, ubj)
Definition mod_tides.F:138
subroutine, public set_iodecomp
Definition set_pio.F:762
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3