1 | MODULE oce |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE oce *** |
---|
4 | !! Ocean : dynamics and active tracers defined in memory |
---|
5 | !!====================================================================== |
---|
6 | !! History : 1.0 ! 2002-11 (G. Madec) F90: Free form and module |
---|
7 | !! 3.1 ! 2009-02 (G. Madec, M. Leclair) pure z* coordinate |
---|
8 | !! 3.3 ! 2010-09 (C. Ethe) TRA-TRC merge: add ts, gtsu, gtsv 4D arrays |
---|
9 | !! 3.7 ! 2014-01 (G. Madec) suppression of curl and before hdiv from in-core memory |
---|
10 | !! 4.1 ! 2019-08 (A. Coward, D. Storkey) rename prognostic variables in preparation for new time scheme |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | USE par_oce ! ocean parameters |
---|
13 | USE lib_mpp ! MPP library |
---|
14 | |
---|
15 | IMPLICIT NONE |
---|
16 | PRIVATE |
---|
17 | |
---|
18 | PUBLIC oce_alloc ! routine called by nemo_init in nemogcm.F90 |
---|
19 | |
---|
20 | !! dynamics and tracer fields |
---|
21 | !! -------------------------- |
---|
22 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: uu , vv !: horizontal velocities [m/s] |
---|
23 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ww !: vertical velocity [m/s] |
---|
24 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: wi !: vertical vel. (adaptive-implicit) [m/s] |
---|
25 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv !: horizontal divergence [s-1] |
---|
26 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) :: ts !: 4D T-S fields [Celsius,psu] |
---|
27 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: rab_b, rab_n !: thermal/haline expansion coef. [Celsius-1,psu-1] |
---|
28 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: rn2b , rn2 !: brunt-vaisala frequency**2 [s-2] |
---|
29 | ! |
---|
30 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: rhd !: in situ density anomalie rhd=(rho-rho0)/rho0 [no units] |
---|
31 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: rhop !: potential volumic mass [kg/m3] |
---|
32 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: Cu_adv !: vertical Courant number (adaptive-implicit) |
---|
33 | |
---|
34 | !! free surface |
---|
35 | !! ------------ |
---|
36 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ssh, uu_b, vv_b !: SSH [m] and barotropic velocities [m/s] |
---|
37 | |
---|
38 | !! Arrays at barotropic time step: ! befbefore! before ! now ! after ! |
---|
39 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ubb_e , ub_e , un_e , ua_e !: u-external velocity |
---|
40 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: vbb_e , vb_e , vn_e , va_e !: v-external velocity |
---|
41 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: sshbb_e, sshb_e, sshn_e, ssha_e !: external ssh |
---|
42 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hu_e !: external u-depth |
---|
43 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hv_e !: external v-depth |
---|
44 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hur_e !: inverse of u-depth |
---|
45 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hvr_e !: inverse of v-depth |
---|
46 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ub2_b , vb2_b !: Half step fluxes (ln_bt_fw=T) |
---|
47 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: un_bf , vn_bf !: Asselin filtered half step fluxes (ln_bt_fw=T) |
---|
48 | #if defined key_agrif |
---|
49 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ub2_i_b, vb2_i_b !: Half step time integrated fluxes |
---|
50 | #endif |
---|
51 | ! |
---|
52 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: spgu, spgv !: horizontal surface pressure gradient |
---|
53 | |
---|
54 | !! interpolated gradient (only used in zps case) |
---|
55 | !! --------------------- |
---|
56 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: gtsu, gtsv !: horizontal gradient of T, S bottom u-point |
---|
57 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: gru , grv !: horizontal gradient of rd at bottom u-point |
---|
58 | |
---|
59 | !! (ISF) interpolated gradient (only used for ice shelf case) |
---|
60 | !! --------------------- |
---|
61 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: gtui, gtvi !: horizontal gradient of T, S and rd at top u-point |
---|
62 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: grui, grvi !: horizontal gradient of T, S and rd at top v-point |
---|
63 | !! (ISF) ice load |
---|
64 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: riceload |
---|
65 | |
---|
66 | !! Energy budget of the leads (open water embedded in sea ice) |
---|
67 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: fraqsr_1lev !: fraction of solar net radiation absorbed in the first ocean level [-] |
---|
68 | INTEGER, PUBLIC, DIMENSION(2) :: noce_array !: unused array but seems to be needed to prevent agrif from creating an empty module |
---|
69 | |
---|
70 | !!---------------------------------------------------------------------- |
---|
71 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
72 | !! $Id$ |
---|
73 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
74 | !!---------------------------------------------------------------------- |
---|
75 | CONTAINS |
---|
76 | |
---|
77 | INTEGER FUNCTION oce_alloc() |
---|
78 | !!---------------------------------------------------------------------- |
---|
79 | !! *** FUNCTION oce_alloc *** |
---|
80 | !!---------------------------------------------------------------------- |
---|
81 | INTEGER :: ierr(6) |
---|
82 | !!---------------------------------------------------------------------- |
---|
83 | ! |
---|
84 | ierr(:) = 0 |
---|
85 | ALLOCATE( uu (jpi,jpj,jpk,jpt) , vv (jpi,jpj,jpk,jpt) , & |
---|
86 | & ww (jpi,jpj,jpk) , hdiv(jpi,jpj,jpk) , & |
---|
87 | & ts (jpi,jpj,jpk,jpts,jpt), & |
---|
88 | & rab_b(jpi,jpj,jpk,jpts) , rab_n(jpi,jpj,jpk,jpts) , & |
---|
89 | & rn2b (jpi,jpj,jpk) , rn2 (jpi,jpj,jpk) , & |
---|
90 | & rhd (jpi,jpj,jpk) , rhop (jpi,jpj,jpk) , STAT=ierr(1) ) |
---|
91 | ! |
---|
92 | ALLOCATE( ssh(jpi,jpj,jpt) , uu_b(jpi,jpj,jpt) , vv_b(jpi,jpj,jpt) , & |
---|
93 | & spgu (jpi,jpj) , spgv(jpi,jpj) , & |
---|
94 | & gtsu(jpi,jpj,jpts), gtsv(jpi,jpj,jpts) , & |
---|
95 | & gru(jpi,jpj) , grv(jpi,jpj) , & |
---|
96 | & gtui(jpi,jpj,jpts), gtvi(jpi,jpj,jpts) , & |
---|
97 | & grui(jpi,jpj) , grvi(jpi,jpj) , & |
---|
98 | & riceload(jpi,jpj) , STAT=ierr(2) ) |
---|
99 | ! |
---|
100 | ALLOCATE( fraqsr_1lev(jpi,jpj) , STAT=ierr(3) ) |
---|
101 | ! |
---|
102 | ALLOCATE( ssha_e(jpi,jpj), sshn_e(jpi,jpj), sshb_e(jpi,jpj), sshbb_e(jpi,jpj), & |
---|
103 | & ua_e(jpi,jpj), un_e(jpi,jpj), ub_e(jpi,jpj), ubb_e(jpi,jpj), & |
---|
104 | & va_e(jpi,jpj), vn_e(jpi,jpj), vb_e(jpi,jpj), vbb_e(jpi,jpj), & |
---|
105 | & hu_e(jpi,jpj), hur_e(jpi,jpj), hv_e(jpi,jpj), hvr_e(jpi,jpj), STAT=ierr(4) ) |
---|
106 | ! |
---|
107 | ALLOCATE( ub2_b(jpi,jpj), vb2_b(jpi,jpj), un_bf(jpi,jpj), vn_bf(jpi,jpj) , STAT=ierr(6) ) |
---|
108 | #if defined key_agrif |
---|
109 | ALLOCATE( ub2_i_b(jpi,jpj), vb2_i_b(jpi,jpj) , STAT=ierr(6) ) |
---|
110 | #endif |
---|
111 | ! |
---|
112 | oce_alloc = MAXVAL( ierr ) |
---|
113 | IF( oce_alloc /= 0 ) CALL ctl_stop( 'STOP', 'oce_alloc: failed to allocate arrays' ) |
---|
114 | ! |
---|
115 | END FUNCTION oce_alloc |
---|
116 | |
---|
117 | !!====================================================================== |
---|
118 | END MODULE oce |
---|