ROMS
Loading...
Searching...
No Matches
mod_arrays Module Reference

Functions/Subroutines

subroutine, public roms_allocate_arrays (allocate_vars)
 
subroutine, public roms_deallocate_arrays
 
subroutine, public roms_initialize_arrays
 

Variables

logical lallocateclima = .TRUE.
 

Function/Subroutine Documentation

◆ roms_allocate_arrays()

subroutine, public mod_arrays::roms_allocate_arrays ( logical, intent(in) allocate_vars)

Definition at line 113 of file mod_arrays.F.

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
220 CALL allocate_nesting
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
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

References mod_average::allocate_average(), mod_bbl::allocate_bbl(), mod_boundary::allocate_boundary(), mod_clima::allocate_clima(), mod_coupling::allocate_coupling(), mod_diags::allocate_diags(), mod_forces::allocate_forces(), mod_grid::allocate_grid(), mod_ice::allocate_ice(), mod_mixing::allocate_mixing(), mod_nesting::allocate_nesting(), mod_ocean::allocate_ocean(), mod_sedbed::allocate_sedbed(), mod_sources::allocate_sources(), mod_tides::allocate_tides(), mod_param::bounds, mod_scalars::extractflag, mod_param::inlm, lallocateclima, mod_scalars::lallocatedmemory, mod_scalars::lclimatology, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, mod_parallel::myrank, mod_param::ngrids, set_pio_mod::set_iodecomp(), wclock_off(), and wclock_on().

Referenced by roms_kernel_mod::roms_initialize(), and roms_kernel_mod::roms_initializep1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ roms_deallocate_arrays()

subroutine, public mod_arrays::roms_deallocate_arrays

Definition at line 257 of file mod_arrays.F.

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
327 CALL deallocate_nesting
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!
336 CALL deallocate_fourdvar
337#endif
338!
339! Deallocate I/O derived-type structures.
340!
341 CALL deallocate_iounits
342!
343! Deallocate main configuration dimensions and associated parameters.
344!
346!
347 RETURN
subroutine, public deallocate_param
Definition mod_param.F:892
integer ngrids
Definition mod_param.F:113

References mod_average::deallocate_average(), mod_bbl::deallocate_bbl(), mod_boundary::deallocate_boundary(), mod_clima::deallocate_clima(), mod_coupling::deallocate_coupling(), mod_diags::deallocate_diags(), mod_forces::deallocate_forces(), mod_fourdvar::deallocate_fourdvar(), mod_grid::deallocate_grid(), mod_ice::deallocate_ice(), mod_iounits::deallocate_iounits(), mod_mixing::deallocate_mixing(), mod_nesting::deallocate_nesting(), mod_ocean::deallocate_ocean(), mod_param::deallocate_param(), mod_sedbed::deallocate_sedbed(), mod_sources::deallocate_sources(), mod_tides::deallocate_tides(), lallocateclima, mod_scalars::lclimatology, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, and mod_param::ngrids.

Referenced by myroms().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ roms_initialize_arrays()

subroutine, public mod_arrays::roms_initialize_arrays

Definition at line 350 of file mod_arrays.F.

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
389 CALL initialize_nesting
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!
439 CALL initialize_fourdvar
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

References mod_parallel::first_tile, mod_average::initialize_average(), mod_bbl::initialize_bbl(), mod_boundary::initialize_boundary(), mod_clima::initialize_clima(), mod_coupling::initialize_coupling(), mod_diags::initialize_diags(), mod_forces::initialize_forces(), mod_fourdvar::initialize_fourdvar(), mod_grid::initialize_grid(), mod_ice::initialize_ice(), mod_mixing::initialize_mixing(), mod_nesting::initialize_nesting(), mod_ocean::initialize_ocean(), mod_sedbed::initialize_sedbed(), mod_tides::initialize_tides(), lallocateclima, mod_parallel::last_tile, mod_scalars::lclimatology, wclock_off(), and wclock_on().

Referenced by roms_kernel_mod::roms_initialize(), and roms_kernel_mod::roms_initializep1().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ lallocateclima

logical mod_arrays::lallocateclima = .TRUE.

Definition at line 106 of file mod_arrays.F.

106 logical :: LallocateClima = .true.

Referenced by roms_allocate_arrays(), roms_deallocate_arrays(), and roms_initialize_arrays().