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

Data Types

type  t_clima
 

Functions/Subroutines

subroutine, public allocate_clima (ng, lbi, ubi, lbj, ubj)
 
subroutine, public deallocate_clima (ng)
 
subroutine, public initialize_clima (ng, tile)
 

Variables

type(t_clima), dimension(:), allocatable clima
 

Function/Subroutine Documentation

◆ allocate_clima()

subroutine, public mod_clima::allocate_clima ( integer, intent(in) ng,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )

Definition at line 157 of file mod_clima.F.

158!
159!=======================================================================
160! !
161! This routine allocates all variables in the module for all nested !
162! grids. !
163! !
164!=======================================================================
165!
166 USE mod_param
167 USE mod_scalars
168!
169! Imported variable declarations.
170!
171 integer, intent(in) :: ng, LBi, UBi, LBj, UBj
172!
173! Local variable declarations.
174!
175 real(r8) :: size2d
176!
177!-----------------------------------------------------------------------
178! Allocate module variables.
179!-----------------------------------------------------------------------
180!
181 IF (ng.eq.1) allocate ( clima(ngrids) )
182!
183! Set horizontal array size.
184!
185 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
186!
187! Climatology/Nudging arrays.
188!
189 IF (lsshclm(ng)) THEN
190 allocate ( clima(ng) % ssh(lbi:ubi,lbj:ubj) )
191 dmem(ng)=dmem(ng)+size2d
192
193#ifndef ANA_SSH
194 allocate ( clima(ng) % sshG(lbi:ubi,lbj:ubj,2) )
195 dmem(ng)=dmem(ng)+2.0_r8*size2d
196#endif
197 END IF
198!
199 IF (lm2clm(ng)) THEN
200 allocate ( clima(ng) % ubarclm(lbi:ubi,lbj:ubj) )
201 dmem(ng)=dmem(ng)+size2d
202
203 allocate ( clima(ng) % vbarclm(lbi:ubi,lbj:ubj) )
204 dmem(ng)=dmem(ng)+size2d
205
206#ifndef ANA_M2CLIMA
207 allocate ( clima(ng) % ubarclmG(lbi:ubi,lbj:ubj,2) )
208 dmem(ng)=dmem(ng)+2.0_r8*size2d
209
210 allocate ( clima(ng) % vbarclmG(lbi:ubi,lbj:ubj,2) )
211 dmem(ng)=dmem(ng)+2.0_r8*size2d
212#endif
213 END IF
214
215#ifdef SOLVE3D
216!
217 IF (lm3clm(ng)) THEN
218 allocate ( clima(ng) % uclm(lbi:ubi,lbj:ubj,n(ng)) )
219 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
220
221 allocate ( clima(ng) % vclm(lbi:ubi,lbj:ubj,n(ng)) )
222 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
223
224# ifndef ANA_M3CLIMA
225 allocate ( clima(ng) % uclmG(lbi:ubi,lbj:ubj,n(ng),2) )
226 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng),r8)*size2d
227
228 allocate ( clima(ng) % vclmG(lbi:ubi,lbj:ubj,n(ng),2) )
229 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng),r8)*size2d
230# endif
231 END IF
232!
233 IF (any(ltracerclm(:,ng)).or.any(lnudgetclm(:,ng))) THEN
234 allocate ( clima(ng) % tclm(lbi:ubi,lbj:ubj,n(ng),ntclm(ng)) )
235 dmem(ng)=dmem(ng)+real(n(ng)*ntclm(ng),r8)*size2d
236
237# ifndef ANA_TCLIMA
238 allocate ( clima(ng) % tclmG(lbi:ubi,lbj:ubj,n(ng),2, &
239 & ntclm(ng)) )
240 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)*ntclm(ng),r8)*size2d
241# endif
242 END IF
243#endif
244!
245! Nudging coefficient arrays.
246!
247 IF (lnudgem2clm(ng)) THEN
248 allocate ( clima(ng) % M2nudgcof(lbi:ubi,lbj:ubj) )
249 dmem(ng)=dmem(ng)+size2d
250 END IF
251
252#ifdef SOLVE3D
253 IF (lnudgem3clm(ng)) THEN
254 allocate ( clima(ng) % M3nudgcof(lbi:ubi,lbj:ubj,n(ng)) )
255 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
256 END IF
257
258 IF (any(lnudgetclm(:,ng))) THEN
259 allocate ( clima(ng) % Tnudgcof(lbi:ubi,lbj:ubj,n(ng), &
260 & ntclm(ng)) )
261 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)*ntclm(ng),r8)*size2d
262 END IF
263#endif
264
265#if defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
266 defined opt_observations || defined sensitivity_4dvar || \
267 defined so_semi
268!
269! Adjoint-based algorithms arrays.
270!
271 allocate ( clima(ng) % zeta_ads(lbi:ubi,lbj:ubj) )
272 dmem(ng)=dmem(ng)+size2d
273
274 allocate ( clima(ng) % zeta_adsG(lbi:ubi,lbj:ubj,2) )
275 dmem(ng)=dmem(ng)+2.0_r8*size2d
276
277 allocate ( clima(ng) % ubar_ads(lbi:ubi,lbj:ubj) )
278 dmem(ng)=dmem(ng)+size2d
279
280 allocate ( clima(ng) % vbar_ads(lbi:ubi,lbj:ubj) )
281 dmem(ng)=dmem(ng)+size2d
282
283 allocate ( clima(ng) % ubar_adsG(lbi:ubi,lbj:ubj,2) )
284 dmem(ng)=dmem(ng)+2.0_r8*size2d
285
286 allocate ( clima(ng) % vbar_adsG(lbi:ubi,lbj:ubj,2) )
287 dmem(ng)=dmem(ng)+2.0_r8*size2d
288
289# ifdef SOLVE3D
290 allocate ( clima(ng) % u_ads(lbi:ubi,lbj:ubj,n(ng)) )
291 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
292
293 allocate ( clima(ng) % v_ads(lbi:ubi,lbj:ubj,n(ng)) )
294 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
295
296 allocate ( clima(ng) % u_adsG(lbi:ubi,lbj:ubj,n(ng),2) )
297 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng),r8)*size2d
298
299 allocate ( clima(ng) % v_adsG(lbi:ubi,lbj:ubj,n(ng),2) )
300 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng),r8)*size2d
301
302 allocate ( clima(ng) % wvel_ads(lbi:ubi,lbj:ubj,0:n(ng)) )
303 dmem(ng)=dmem(ng)+real(n(ng)+1,r8)*size2d
304
305 allocate ( clima(ng) % wvel_adsG(lbi:ubi,lbj:ubj,0:n(ng),2) )
306 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)+1,r8)*size2d
307
308 allocate ( clima(ng) % t_ads(lbi:ubi,lbj:ubj,n(ng),nt(ng)) )
309 dmem(ng)=dmem(ng)+real(n(ng)*nt(ng),r8)*size2d
310
311 allocate ( clima(ng) % t_adsG(lbi:ubi,lbj:ubj,n(ng),2,nt(ng)) )
312 dmem(ng)=dmem(ng)+2.0_r8*real(n(ng)*nt(ng),r8)*size2d
313# endif
314#endif
315!
316 RETURN
integer, dimension(:), allocatable ntclm
Definition mod_param.F:494
integer, dimension(:), allocatable n
Definition mod_param.F:479
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
integer ngrids
Definition mod_param.F:113
integer, dimension(:), allocatable nt
Definition mod_param.F:489
logical, dimension(:), allocatable lnudgem2clm
logical, dimension(:), allocatable lm3clm
logical, dimension(:), allocatable lsshclm
logical, dimension(:), allocatable lnudgem3clm
logical, dimension(:), allocatable lm2clm
logical, dimension(:,:), allocatable ltracerclm
logical, dimension(:,:), allocatable lnudgetclm

References clima, mod_param::dmem, mod_scalars::lm2clm, mod_scalars::lm3clm, mod_scalars::lnudgem2clm, mod_scalars::lnudgem3clm, mod_scalars::lnudgetclm, mod_scalars::lsshclm, mod_scalars::ltracerclm, mod_param::n, mod_param::ngrids, mod_param::nt, mod_param::ntclm, and mod_kinds::r8.

Referenced by mod_arrays::roms_allocate_arrays().

Here is the caller graph for this function:

◆ deallocate_clima()

subroutine, public mod_clima::deallocate_clima ( integer, intent(in) ng)

Definition at line 319 of file mod_clima.F.

320!
321!=======================================================================
322! !
323! This routine deallocates all variables in the module for all nested !
324! grids. !
325! !
326!=======================================================================
327!
328 USE mod_param
329 USE mod_scalars
330
331#ifdef SUBOBJECT_DEALLOCATION
332!
333 USE destroy_mod, ONLY : destroy
334#endif
335!
336! Imported variable declarations.
337!
338 integer, intent(in) :: ng
339!
340! Local variable declarations.
341!
342 character (len=*), parameter :: MyFile = &
343 & __FILE__//", deallocate_clima"
344
345# ifdef SUBOBJECT_DEALLOCATION
346!
347!-----------------------------------------------------------------------
348! Deallocate each variable in the derived-type T_CLIMA structure
349! separately.
350!-----------------------------------------------------------------------
351!
352! Climatology/Nudging arrays.
353!
354 IF (lsshclm(ng)) THEN
355 IF (.not.destroy(ng, clima(ng)%ssh, myfile, &
356 & __line__, 'CLIMA(ng)%ssh')) RETURN
357
358# ifndef ANA_SSH
359 IF (.not.destroy(ng, clima(ng)%sshG, myfile, &
360 & __line__, 'CLIMA(ng)%sshG')) RETURN
361# endif
362 END IF
363!
364 IF (lm2clm(ng)) THEN
365 IF (.not.destroy(ng, clima(ng)%ubarclm, myfile, &
366 & __line__, 'CLIMA(ng)%ubarclm')) RETURN
367
368 IF (.not.destroy(ng, clima(ng)%vbarclm, myfile, &
369 & __line__, 'CLIMA(ng)%vbarclm')) RETURN
370
371# ifndef ANA_M2CLIMA
372 IF (.not.destroy(ng, clima(ng)%ubarclmG, myfile, &
373 & __line__, 'CLIMA(ng)%ubarclmG')) RETURN
374
375 IF (.not.destroy(ng, clima(ng)%vbarclmG, myfile, &
376 & __line__, 'CLIMA(ng)%vbarclmG')) RETURN
377# endif
378 END IF
379
380# ifdef SOLVE3D
381!
382 IF (lm3clm(ng)) THEN
383 IF (.not.destroy(ng, clima(ng)%uclm, myfile, &
384 & __line__, 'CLIMA(ng)%uclm')) RETURN
385
386 IF (.not.destroy(ng, clima(ng)%vclm, myfile, &
387 & __line__, 'CLIMA(ng)%vclm')) RETURN
388
389# ifndef ANA_M3CLIMA
390 IF (.not.destroy(ng, clima(ng)%uclmG, myfile, &
391 & __line__, 'CLIMA(ng)%uclmG')) RETURN
392
393 IF (.not.destroy(ng, clima(ng)%vclmG, myfile, &
394 & __line__, 'CLIMA(ng)%vclmG')) RETURN
395# endif
396 END IF
397!
398 IF (any(ltracerclm(:,ng)).or.any(lnudgetclm(:,ng))) THEN
399 IF (.not.destroy(ng, clima(ng)%tclm, myfile, &
400 & __line__, 'CLIMA(ng)%tclm')) RETURN
401
402# ifndef ANA_TCLIMA
403 IF (.not.destroy(ng, clima(ng)%tclmG, myfile, &
404 & __line__, 'CLIMA(ng)%tclmG')) RETURN
405# endif
406 END IF
407# endif
408!
409! Nudging coefficient arrays.
410!
411 IF (lnudgem2clm(ng)) THEN
412 IF (.not.destroy(ng, clima(ng)%M2nudgcof, myfile, &
413 & __line__, 'CLIMA(ng)%M2nudgcof')) RETURN
414 END IF
415
416# ifdef SOLVE3D
417 IF (lnudgem3clm(ng)) THEN
418 IF (.not.destroy(ng, clima(ng)%M3nudgcof, myfile, &
419 & __line__, 'CLIMA(ng)%M3nudgcof')) RETURN
420 END IF
421
422 IF (any(lnudgetclm(:,ng))) THEN
423 IF (.not.destroy(ng, clima(ng)%Tnudgcof, myfile, &
424 & __line__, 'CLIMA(ng)%Tnudgcof')) RETURN
425 END IF
426# endif
427
428# if defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
429 defined opt_observations || defined sensitivity_4dvar || \
430 defined so_semi
431!
432! Adjoint-based algorithms arrays.
433!
434 IF (.not.destroy(ng, clima(ng)%zeta_ads, myfile, &
435 & __line__, 'CLIMA(ng)%zeta_ads')) RETURN
436
437 IF (.not.destroy(ng, clima(ng)%zeta_adsG, myfile, &
438 & __line__, 'CLIMA(ng)%zeta_adsG')) RETURN
439
440 IF (.not.destroy(ng, clima(ng)%ubar_ads, myfile, &
441 & __line__, 'CLIMA(ng)%ubar_ads')) RETURN
442
443 IF (.not.destroy(ng, clima(ng)%vbar_ads, myfile, &
444 & __line__, 'CLIMA(ng)%vbar_ads')) RETURN
445
446 IF (.not.destroy(ng, clima(ng)%ubar_adsG, myfile, &
447 & __line__, 'CLIMA(ng)%ubar_adsG')) RETURN
448
449 IF (.not.destroy(ng, clima(ng)%vbar_adsG, myfile, &
450 & __line__, 'CLIMA(ng)%vbar_adsG')) RETURN
451
452# ifdef SOLVE3D
453 IF (.not.destroy(ng, clima(ng)%u_ads, myfile, &
454 & __line__, 'CLIMA(ng)%u_ads')) RETURN
455
456 IF (.not.destroy(ng, clima(ng)%v_ads, myfile, &
457 & __line__, 'CLIMA(ng)%v_ads')) RETURN
458
459 IF (.not.destroy(ng, clima(ng)%u_adsG, myfile, &
460 & __line__, 'CLIMA(ng)%u_adsG')) RETURN
461
462 IF (.not.destroy(ng, clima(ng)%v_adsG, myfile, &
463 & __line__, 'CLIMA(ng)%v_adsG')) RETURN
464
465 IF (.not.destroy(ng, clima(ng)%wvel_ads, myfile, &
466 & __line__, 'CLIMA(ng)%wvel_ads')) RETURN
467
468 IF (.not.destroy(ng, clima(ng)%wvel_adsG, myfile, &
469 & __line__, 'CLIMA(ng)%wvel_adsG')) RETURN
470
471 IF (.not.destroy(ng, clima(ng)%t_ads, myfile, &
472 & __line__, 'CLIMA(ng)%t_ads')) RETURN
473
474 IF (.not.destroy(ng, clima(ng)%t_adsG, myfile, &
475 & __line__, 'CLIMA(ng)%t_adsG')) RETURN
476# endif
477# endif
478#endif
479!
480!-----------------------------------------------------------------------
481! Deallocate derived-type CLIMA structure.
482!-----------------------------------------------------------------------
483!
484 IF (ng.eq.ngrids) THEN
485 IF (allocated(clima)) deallocate ( clima )
486 END IF
487!
488 RETURN

References clima, mod_scalars::lm2clm, mod_scalars::lm3clm, mod_scalars::lnudgem2clm, mod_scalars::lnudgem3clm, mod_scalars::lnudgetclm, mod_scalars::lsshclm, mod_scalars::ltracerclm, and mod_param::ngrids.

Referenced by mod_arrays::roms_deallocate_arrays().

Here is the caller graph for this function:

◆ initialize_clima()

subroutine, public mod_clima::initialize_clima ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 491 of file mod_clima.F.

492!
493!=======================================================================
494! !
495! This routine initialize all variables in the module using first !
496! touch distribution policy. In shared-memory configuration, this !
497! operation actually performs propagation of the "shared arrays" !
498! across the cluster, unless another policy is specified to !
499! override the default. !
500! !
501!=======================================================================
502!
503 USE mod_param
504 USE mod_scalars
505!
506! Imported variable declarations.
507!
508 integer, intent(in) :: ng, tile
509!
510! Local variable declarations.
511!
512 integer :: Imin, Imax, Jmin, Jmax
513 integer :: i, j
514#ifdef SOLVE3D
515 integer :: itrc, k
516#endif
517
518 real(r8), parameter :: IniVal = 0.0_r8
519
520#include "set_bounds.h"
521!
522! Set array initialization range.
523!
524#ifdef DISTRIBUTE
525 imin=bounds(ng)%LBi(tile)
526 imax=bounds(ng)%UBi(tile)
527 jmin=bounds(ng)%LBj(tile)
528 jmax=bounds(ng)%UBj(tile)
529#else
530 IF (domain(ng)%Western_Edge(tile)) THEN
531 imin=bounds(ng)%LBi(tile)
532 ELSE
533 imin=istr
534 END IF
535 IF (domain(ng)%Eastern_Edge(tile)) THEN
536 imax=bounds(ng)%UBi(tile)
537 ELSE
538 imax=iend
539 END IF
540 IF (domain(ng)%Southern_Edge(tile)) THEN
541 jmin=bounds(ng)%LBj(tile)
542 ELSE
543 jmin=jstr
544 END IF
545 IF (domain(ng)%Northern_Edge(tile)) THEN
546 jmax=bounds(ng)%UBj(tile)
547 ELSE
548 jmax=jend
549 END IF
550#endif
551!
552!-----------------------------------------------------------------------
553! Initialize module variables.
554!-----------------------------------------------------------------------
555!
556! Climatology/Nudging arrays.
557!
558 IF (lsshclm(ng)) THEN
559 DO j=jmin,jmax
560 DO i=imin,imax
561 clima(ng) % ssh(i,j) = inival
562#ifndef ANA_SSH
563 clima(ng) % sshG(i,j,1) = inival
564 clima(ng) % sshG(i,j,2) = inival
565#endif
566 END DO
567 END DO
568 END IF
569!
570 IF (lm2clm(ng)) THEN
571 DO j=jmin,jmax
572 DO i=imin,imax
573 clima(ng) % ubarclm(i,j) = inival
574 clima(ng) % vbarclm(i,j) = inival
575#ifndef ANA_M2CLIMA
576 clima(ng) % ubarclmG(i,j,1) = inival
577 clima(ng) % ubarclmG(i,j,2) = inival
578 clima(ng) % vbarclmG(i,j,1) = inival
579 clima(ng) % vbarclmG(i,j,2) = inival
580#endif
581 END DO
582 END DO
583 END IF
584
585#ifdef SOLVE3D
586!
587 IF (lm3clm(ng)) THEN
588 DO k=1,n(ng)
589 DO j=jmin,jmax
590 DO i=imin,imax
591 clima(ng) % uclm(i,j,k) = inival
592 clima(ng) % vclm(i,j,k) = inival
593# ifndef ANA_M3CLIMA
594 clima(ng) % uclmG(i,j,k,1) = inival
595 clima(ng) % uclmG(i,j,k,2) = inival
596 clima(ng) % vclmG(i,j,k,1) = inival
597 clima(ng) % vclmG(i,j,k,2) = inival
598# endif
599 END DO
600 END DO
601 END DO
602 END IF
603!
604 IF (any(ltracerclm(:,ng)).or.any(lnudgetclm(:,ng))) THEN
605 DO itrc=1,ntclm(ng)
606 DO k=1,n(ng)
607 DO j=jmin,jmax
608 DO i=imin,imax
609 clima(ng) % tclm(i,j,k,itrc) = inival
610# ifndef ANA_TCLIMA
611 clima(ng) % tclmG(i,j,k,1,itrc) = inival
612 clima(ng) % tclmG(i,j,k,2,itrc) = inival
613# endif
614 END DO
615 END DO
616 END DO
617 END DO
618 END IF
619#endif
620!
621! Nudging coefficient arrays.
622!
623 IF (lnudgem2clm(ng)) THEN
624 DO j=jmin,jmax
625 DO i=imin,imax
626 clima(ng) % M2nudgcof(i,j) = inival
627 END DO
628 END DO
629 END IF
630
631#ifdef SOLVE3D
632!
633 IF (lnudgem3clm(ng)) THEN
634 DO k=1,n(ng)
635 DO j=jmin,jmax
636 DO i=imin,imax
637 clima(ng) % M3nudgcof(i,j,k) = inival
638 END DO
639 END DO
640 END DO
641 END IF
642!
643 IF (any(lnudgetclm(:,ng))) THEN
644 DO itrc=1,ntclm(ng)
645 DO k=1,n(ng)
646 DO j=jmin,jmax
647 DO i=imin,imax
648 clima(ng) % Tnudgcof(i,j,k,itrc) = inival
649 END DO
650 END DO
651 END DO
652 END DO
653 END IF
654#endif
655
656#if defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
657 defined opt_observations || defined sensitivity_4dvar || \
658 defined so_semi
659!
660! Adjoint-based algorithms arrays.
661!
662 DO j=jmin,jmax
663 DO i=imin,imax
664 clima(ng) % zeta_ads(i,j) = inival
665 clima(ng) % zeta_adsG(i,j,1) = inival
666 clima(ng) % zeta_adsG(i,j,2) = inival
667!
668 clima(ng) % ubar_ads(i,j) = inival
669 clima(ng) % vbar_ads(i,j) = inival
670 clima(ng) % ubar_adsG(i,j,1) = inival
671 clima(ng) % ubar_adsG(i,j,2) = inival
672 clima(ng) % vbar_adsG(i,j,1) = inival
673 clima(ng) % vbar_adsG(i,j,2) = inival
674 END DO
675
676# ifdef SOLVE3D
677!
678 DO k=1,n(ng)
679 DO i=imin,imax
680 clima(ng) % u_ads(i,j,k) = inival
681 clima(ng) % v_ads(i,j,k) = inival
682 clima(ng) % u_adsG(i,j,k,1) = inival
683 clima(ng) % u_adsG(i,j,k,2) = inival
684 clima(ng) % v_adsG(i,j,k,1) = inival
685 clima(ng) % v_adsG(i,j,k,2) = inival
686 END DO
687 END DO
688!
689 DO k=0,n(ng)
690 DO i=imin,imax
691 clima(ng) % wvel_ads(i,j,k) = inival
692 clima(ng) % wvel_adsG(i,j,k,1) = inival
693 clima(ng) % wvel_adsG(i,j,k,2) = inival
694 END DO
695 END DO
696!
697 DO itrc=1,nt(ng)
698 DO k=1,n(ng)
699 DO i=imin,imax
700 clima(ng) % t_ads(i,j,k,itrc) = inival
701 clima(ng) % t_adsG(i,j,k,1,itrc) = inival
702 clima(ng) % t_adsG(i,j,k,2,itrc) = inival
703 END DO
704 END DO
705 END DO
706# endif
707 END DO
708#endif
709!
710 RETURN
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329

References mod_param::bounds, clima, mod_param::domain, mod_scalars::lm2clm, mod_scalars::lm3clm, mod_scalars::lnudgem2clm, mod_scalars::lnudgem3clm, mod_scalars::lnudgetclm, mod_scalars::lsshclm, mod_scalars::ltracerclm, mod_param::n, mod_param::nt, and mod_param::ntclm.

Referenced by mod_arrays::roms_initialize_arrays().

Here is the caller graph for this function:

Variable Documentation

◆ clima

type (t_clima), dimension(:), allocatable mod_clima::clima

Definition at line 153 of file mod_clima.F.

153 TYPE (T_CLIMA), allocatable :: CLIMA(:)

Referenced by ad_nesting_mod::ad_correct_tracer_tile(), ad_get_data(), ad_rhs3d_mod::ad_rhs3d_tile(), ad_set_data_tile(), ad_step2d_mod::ad_step2d_tile(), ad_step3d_t_mod::ad_step3d_t_tile(), ad_t3dbc_mod::ad_t3dbc_tile(), ad_t3dmix2_mod::ad_t3dmix2(), ad_t3dmix4_mod::ad_t3dmix4(), ad_u2dbc_mod::ad_u2dbc_tile(), ad_u3dbc_mod::ad_u3dbc_tile(), ad_v2dbc_mod::ad_v2dbc_tile(), ad_v3dbc_mod::ad_v3dbc_tile(), adsen_force_mod::adsen_force(), adsen_initial_mod::adsen_initial(), allocate_clima(), analytical_mod::ana_m2clima_tile(), analytical_mod::ana_m3clima_tile(), analytical_mod::ana_nudgcoef_tile(), analytical_mod::ana_ssh_tile(), analytical_mod::ana_tclima_tile(), nesting_mod::correct_tracer_tile(), deallocate_clima(), get_data(), get_nudgcoef_mod::get_nudgcoef_nf90(), get_nudgcoef_mod::get_nudgcoef_pio(), ice_thermo_mod::ice_thermo(), initialize_clima(), rhs3d_mod::rhs3d_tile(), rp_get_data(), rp_rhs3d_mod::rp_rhs3d_tile(), rp_set_data_tile(), rp_step2d_mod::rp_step2d_tile(), rp_step3d_t_mod::rp_step3d_t_tile(), rp_t3dbc_mod::rp_t3dbc_tile(), rp_t3dmix2_mod::rp_t3dmix2(), rp_t3dmix4_mod::rp_t3dmix4(), rp_u2dbc_mod::rp_u2dbc_tile(), rp_u3dbc_mod::rp_u3dbc_tile(), rp_v2dbc_mod::rp_v2dbc_tile(), rp_v3dbc_mod::rp_v3dbc_tile(), set_data_tile(), set_tides_mod::set_tides_tile(), step2d_mod::step2d_tile(), step3d_t_mod::step3d_t_tile(), t3dbc_mod::t3dbc_tile(), t3dmix2_mod::t3dmix2(), t3dmix4_mod::t3dmix4(), tl_nesting_mod::tl_correct_tracer_tile(), tl_get_data(), tl_rhs3d_mod::tl_rhs3d_tile(), tl_set_data_tile(), tl_step2d_mod::tl_step2d_tile(), tl_step3d_t_mod::tl_step3d_t_tile(), tl_t3dbc_mod::tl_t3dbc_tile(), tl_t3dmix2_mod::tl_t3dmix2(), tl_t3dmix4_mod::tl_t3dmix4(), tl_u2dbc_mod::tl_u2dbc_tile(), tl_u3dbc_mod::tl_u3dbc_tile(), tl_v2dbc_mod::tl_v2dbc_tile(), tl_v3dbc_mod::tl_v3dbc_tile(), u2dbc_mod::u2dbc_tile(), u3dbc_mod::u3dbc_tile(), v2dbc_mod::v2dbc_tile(), and v3dbc_mod::v3dbc_tile().