New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
oce_tam.F90 in NEMO/branches/NERC/dev_release-3.4_NEMOTAM_consolidated/NEMOGCM/NEMO/OPATAM_SRC – NEMO

source: NEMO/branches/NERC/dev_release-3.4_NEMOTAM_consolidated/NEMOGCM/NEMO/OPATAM_SRC/oce_tam.F90 @ 13432

Last change on this file since 13432 was 13420, checked in by smueller, 4 years ago

Addition of interior masks as suggested in ticket #1499

  • Property svn:executable set to *
File size: 19.8 KB
Line 
1MODULE oce_tam
2   !!----------------------------------------------------------------------
3   !!    This software is governed by the CeCILL licence (Version 2)
4   !!----------------------------------------------------------------------
5   !!======================================================================
6   !!                       ***  MODULE oce_tam ***
7   !! NEMOVAR Tangent linear and Adjoint Model variables :
8   !!
9   !!    Allocate tangent linear and adjoint fields for the inner loop
10   !!
11   !!======================================================================
12   !! History of the direct module:
13   !!   8.5  !  02-11  (G. Madec)  F90: Free form and module
14   !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
15   !! History of the TAM module:
16   !!   9.0  !  07-07 (K. Mogensen) Initial version
17   !!   9.0  !  08-03 (A. Vidard) Add variables
18   !!   9.0  !  09-03 (A. Weaver) Nemo v3 compatible, merge tl_init/ad_init
19   !!        ! 2011-07 (D. Lea) Add altimeter bias and sea ice
20   !! NEMO 3.4  ! 2012-04 (P.-A. Bouttier) update 3.4
21   !!           ! 2012-09 (A. Vidard) Deallocating and initialising options
22   !!----------------------------------------------------------------------
23   !!   oce_tam_init : Allocate and initialize the TAM fields
24   !!----------------------------------------------------------------------
25   !! * Modules used
26   USE par_oce
27   USE lib_mpp
28   USE dom_oce
29   USE domwri
30
31   IMPLICIT NONE
32
33   !! * Routine accessibility
34   PRIVATE
35
36   !! dynamics and tracer fields                            ! before ! now    ! after  ! the after trends becomes the fields
37   !! --------------------------                            ! fields ! fields ! trends ! only after tra_zdf and dyn_spg
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ub_tl   ,  un_tl    , ua_tl     !: i-horizontal velocity        [m/s]
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   vb_tl   ,  vn_tl    , va_tl     !: j-horizontal velocity        [m/s]
40   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::              wn_tl                !: vertical velocity            [m/s]
41   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rotb_tl ,  rotn_tl              !: relative vorticity           [s-1]
42   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   hdivb_tl,  hdivn_tl             !: horizontal divergence        [s-1]
43   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsb_tl  ,  tsn_tl               !: 4D T-S fields        [Celcius,psu]
44   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rn2b_tl ,  rn2_tl               !: brunt-vaisala frequency**2   [s-2]
45   !
46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::  tsa_tl                           !: 4D T-S trends fields & work array
47   !
48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhd_tl                            !: in situ density anomalie rhd=(rho-rau0)/rau0  [no units]
49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhop_tl                           !: potential volumic mass                           [kg/m3]
50
51   !! free surface                                      !  before  ! now    ! after  !
52   !! ------------                                      !  fields  ! fields ! trends !
53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET ::   sshb_tl   , sshn_tl   , ssha_tl   !: sea surface height at t-point [m]
54   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::   sshu_b_tl , sshu_n_tl , sshu_a_tl !: sea surface height at u-point [m]
55   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::   sshv_b_tl , sshv_n_tl , sshv_a_tl !: sea surface height at u-point [m]
56   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::            sshf_n_tl                !: sea surface height at f-point [m]
57   !
58   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   spgu_tl, spgv_tl                          !: horizontal surface pressure gradient
59
60   !! interpolated gradient (only used in zps case)
61   !! ---------------------
62   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtsu_tl, gtsv_tl   !: horizontal gradient of T, S bottom u-point
63   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gru_tl , grv_tl    !: horizontal gradient of rd at bottom u-point
64   !
65   LOGICAL, SAVE, PRIVATE :: ll_alloctl = .FALSE.
66   !
67   !! dynamics and tracer fields                            ! before ! now    ! after  ! the after trends becomes the fields
68   !! --------------------------                            ! fields ! fields ! trends ! only after tra_zdf and dyn_spg
69   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ub_ad   ,  un_ad    , ua_ad     !: i-horizontal velocity        [m/s]
70   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   vb_ad   ,  vn_ad    , va_ad     !: j-horizontal velocity        [m/s]
71   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::              wn_ad                !: vertical velocity            [m/s]
72   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rotb_ad ,  rotn_ad              !: relative vorticity           [s-1]
73   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   hdivb_ad,  hdivn_ad             !: horizontal divergence        [s-1]
74   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsb_ad  ,  tsn_ad               !: 4D T-S fields        [Celcius,psu]
75   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   rn2b_ad ,  rn2_ad               !: brunt-vaisala frequency**2   [s-2]
76   !
77   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::  tsa_ad                           !: 4D T-S trends fields & work array
78   !
79   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhd_ad                            !: in situ density anomalie rhd=(rho-rau0)/rau0  [no units]
80   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   rhop_ad                           !: potential volumic mass                           [kg/m3]
81
82   !! free surface                                      !  before  ! now    ! after  !
83   !! ------------                                      !  fields  ! fields ! trends !
84   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET ::   sshb_ad   , sshn_ad   , ssha_ad   !: sea surface height at t-point [m]
85   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::   sshu_b_ad , sshu_n_ad , sshu_a_ad !: sea surface height at u-point [m]
86   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::   sshv_b_ad , sshv_n_ad , sshv_a_ad !: sea surface height at u-point [m]
87   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)         ::            sshf_n_ad                !: sea surface height at f-point [m]
88   !
89   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   spgu_ad, spgv_ad                          !: horizontal surface pressure gradient
90
91   !! interpolated gradient (only used in zps case)
92   !! ---------------------
93   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   gtsu_ad, gtsv_ad   !: horizontal gradient of T, S bottom u-point
94   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   gru_ad , grv_ad    !: horizontal gradient of rd at bottom u-point
95
96   LOGICAL, PRIVATE, SAVE :: ll_allocad = .FALSE.
97#if defined key_zdfddm
98!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99!!!!  AW: The declaration/allocation/initialization of these variables
100!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
101!!!!      with NEMO.
102   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE :: &
103      & rrau_tl,  & !: Tangent Linear of heat/salt buoyancy flux ratio
104      & rrau_ad     !: Adjoint of heat/salt buoyancy flux ratio
105!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
106#endif
107   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: &
108      & tmsk_i, &
109      & umsk_i, &
110      & vmsk_i, &
111      & fmsk_i
112
113   !!----------------------------------------------------------------------
114   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
115   !! $Id$
116   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
117   !!----------------------------------------------------------------------
118
119   PUBLIC oce_alloc_tam
120   PUBLIC oce_dealloc_tam
121   PUBLIC oce_tam_init
122
123CONTAINS
124   INTEGER FUNCTION oce_alloc_tam( kmode )
125      !!----------------------------------------------------------------------
126      !!                   ***  FUNCTION oce_alloc  ***
127      !!----------------------------------------------------------------------
128      INTEGER, OPTIONAL :: kmode
129      INTEGER :: ierr(5)
130      INTEGER :: jmode
131      INTEGER :: jk
132      REAL(wp), DIMENSION(jpi,jpj) :: &
133         & z_tmsk_i, &
134         & z_umsk_i, &
135         & z_vmsk_i, &
136         & z_fmsk_i
137
138      !!----------------------------------------------------------------------
139      !
140      IF ( PRESENT(kmode) ) THEN
141         jmode = kmode
142      ELSE
143         jmode = 0
144      END IF
145
146      IF (.NOT. ALLOCATED ( tmsk_i ) ) THEN
147         ALLOCATE ( tmsk_i (jpi, jpj, jpk), &
148            &       umsk_i (jpi, jpj, jpk), &
149            &       vmsk_i (jpi, jpj, jpk), &
150            &       fmsk_i (jpi, jpj, jpk) )
151
152         CALL dom_uniq( z_tmsk_i, 'T', 1 )
153         CALL dom_uniq( z_umsk_i, 'U', 1 )
154         CALL dom_uniq( z_vmsk_i, 'V', 1 )
155         CALL dom_uniq( z_fmsk_i, 'F', 1 )
156
157         DO jk = 1, jpk
158            tmsk_i(:,:,jk) = tmask(:,:,jk) * z_tmsk_i(:,:)
159            umsk_i(:,:,jk) = umask(:,:,jk) * z_umsk_i(:,:)
160            vmsk_i(:,:,jk) = vmask(:,:,jk) * z_vmsk_i(:,:)
161            fmsk_i(:,:,jk) = fmask(:,:,jk) * z_fmsk_i(:,:)
162         END DO
163
164      END IF
165
166      IF ( ( jmode == 0 ) .OR. ( jmode == 1 ) .AND. ( .NOT. ll_alloctl) ) THEN
167      ALLOCATE( ub_tl   (jpi,jpj,jpk)      , un_tl   (jpi,jpj,jpk)      , ua_tl(jpi,jpj,jpk)       ,     &
168         &      vb_tl   (jpi,jpj,jpk)      , vn_tl   (jpi,jpj,jpk)      , va_tl(jpi,jpj,jpk)       ,     &
169         &      wn_tl   (jpi,jpj,jpk)      ,                                                       &
170         &      rotb_tl (jpi,jpj,jpk)      , rotn_tl (jpi,jpj,jpk)      ,                             &
171         &      hdivb_tl(jpi,jpj,jpk)      , hdivn_tl(jpi,jpj,jpk)      ,                             &
172         &      tsb_tl  (jpi,jpj,jpk,jpts) , tsn_tl  (jpi,jpj,jpk,jpts) , tsa_tl(jpi,jpj,jpk,jpts) ,     &
173         &      rn2b_tl (jpi,jpj,jpk)      , rn2_tl  (jpi,jpj,jpk)                              , STAT=ierr(1) )
174         !
175      ALLOCATE(rhd_tl (jpi,jpj,jpk) ,                                         &
176         &     rhop_tl(jpi,jpj,jpk) ,                                         &
177         &     sshb_tl  (jpi,jpj)   , sshn_tl  (jpi,jpj) , ssha_tl  (jpi,jpj) ,     &
178         &     sshu_b_tl(jpi,jpj)   , sshu_n_tl(jpi,jpj) , sshu_a_tl(jpi,jpj) ,     &
179         &     sshv_b_tl(jpi,jpj)   , sshv_n_tl(jpi,jpj) , sshv_a_tl(jpi,jpj) ,     &
180         &                         sshf_n_tl(jpi,jpj) ,                       &
181         &     spgu_tl  (jpi,jpj)   , spgv_tl(jpi,jpj)   ,                       &
182         &     gtsu_tl(jpi,jpj,jpts), gtsv_tl(jpi,jpj,jpts),                     &
183         &     gru_tl(jpi,jpj)      , grv_tl(jpi,jpj)                      , STAT=ierr(2) )
184
185#if defined key_zdfddm
186         ALLOCATE( rrau_tl(jpi,jpj,jpk), STAT=ierr(5) )
187#endif
188         !
189         ll_alloctl = .TRUE.
190      END IF
191      !
192      IF ( ( jmode == 0 ) .OR. ( jmode == 2 ) .AND. ( .NOT. ll_allocad)  ) THEN
193         ALLOCATE( ub_ad   (jpi,jpj,jpk)      , un_ad   (jpi,jpj,jpk)      , ua_ad(jpi,jpj,jpk)       , &
194            &      vb_ad   (jpi,jpj,jpk)      , vn_ad   (jpi,jpj,jpk)      , va_ad(jpi,jpj,jpk)       , &
195            &      wn_ad   (jpi,jpj,jpk)      ,                                                         &
196            &      rotb_ad (jpi,jpj,jpk)      , rotn_ad (jpi,jpj,jpk)      ,                            &
197            &      hdivb_ad(jpi,jpj,jpk)      , hdivn_ad(jpi,jpj,jpk)      ,                            &
198            &      tsb_ad  (jpi,jpj,jpk,jpts) , tsn_ad  (jpi,jpj,jpk,jpts) , tsa_ad(jpi,jpj,jpk,jpts) , &
199            &      rn2b_ad (jpi,jpj,jpk)      , rn2_ad  (jpi,jpj,jpk)                              , STAT=ierr(3) )
200         !
201         ALLOCATE(rhd_ad (jpi,jpj,jpk) ,                                           &
202            &     rhop_ad(jpi,jpj,jpk) ,                                           &
203            &     sshb_ad  (jpi,jpj)   , sshn_ad  (jpi,jpj) , ssha_ad  (jpi,jpj) , &
204            &     sshu_b_ad(jpi,jpj)   , sshu_n_ad(jpi,jpj) , sshu_a_ad(jpi,jpj) , &
205            &     sshv_b_ad(jpi,jpj)   , sshv_n_ad(jpi,jpj) , sshv_a_ad(jpi,jpj) , &
206            &                         sshf_n_ad(jpi,jpj) ,                         &
207            &     spgu_ad  (jpi,jpj)   , spgv_ad(jpi,jpj)   ,                      &
208            &     gtsu_ad(jpi,jpj,jpts), gtsv_ad(jpi,jpj,jpts),                    &
209            &     gru_ad(jpi,jpj)      , grv_ad(jpi,jpj)                    , STAT=ierr(4) )
210
211#if defined key_zdfddm
212         ALLOCATE( rrau_ad(jpi,jpj,jpk), STAT=ierr(5) )
213#endif
214
215         ll_allocad = .TRUE.
216      END IF
217      oce_alloc_tam = MAXVAL( ierr )
218      IF( oce_alloc_tam /= 0 )   CALL ctl_warn('oce_alloc_tam: failed to allocate arrays')
219      !
220   END FUNCTION oce_alloc_tam
221   !
222   INTEGER FUNCTION oce_dealloc_tam( kmode )
223      !!----------------------------------------------------------------------
224      !!                   ***  FUNCTION oce_dealloc  ***
225      !!----------------------------------------------------------------------
226      INTEGER, OPTIONAL :: kmode
227      INTEGER :: jmode
228      INTEGER :: ierr(5)
229      !!----------------------------------------------------------------------
230      !
231      IF ( PRESENT(kmode) ) THEN
232         jmode = kmode
233      ELSE
234         jmode = 0
235      END IF
236
237      IF ( ( jmode == 0 ) .OR. ( jmode == 1 ) .AND. ( ll_alloctl ) ) THEN
238         DEALLOCATE( ub_tl         , un_tl         , ua_tl       , &
239            &        vb_tl         , vn_tl         , va_tl       , &
240            &        wn_tl         ,                               &
241            &        rotb_tl       , rotn_tl       ,               &
242            &        hdivb_tl      , hdivn_tl      ,               &
243            &        tsb_tl        , tsn_tl        , tsa_tl      , &
244            &        rn2b_tl       , rn2_tl                      , STAT=ierr(1) )
245         !
246         DEALLOCATE( rhd_tl     ,                                  &
247            &       rhop_tl     ,                                  &
248            &       sshb_tl     , sshn_tl   , ssha_tl   ,          &
249            &       sshu_b_tl   , sshu_n_tl , sshu_a_tl ,          &
250            &       sshv_b_tl   , sshv_n_tl , sshv_a_tl ,          &
251            &       sshf_n_tl   ,                                  &
252            &       spgu_tl     , spgv_tl   ,                      &
253            &       gtsu_tl     , gtsv_tl   ,                      &
254            &       gru_tl      , grv_tl                         , STAT=ierr(2) )
255
256#if defined key_zdfddm
257         DEALLOCATE( rrau_tl, STAT=ierr(5) )
258#endif
259         !
260         ll_alloctl = .FALSE.
261      END IF
262      !
263      IF ( ( jmode == 0 ) .OR. ( jmode == 2 ) .and. ( ll_allocad ) ) THEN
264         DEALLOCATE( ub_ad       , un_ad         , ua_ad       , &
265            &        vb_ad       , vn_ad         , va_ad       , &
266            &        wn_ad       ,                               &
267            &        rotb_ad     , rotn_ad       ,               &
268            &        hdivb_ad    , hdivn_ad      ,               &
269            &        tsb_ad      , tsn_ad        , tsa_ad      , &
270            &        rn2b_ad     , rn2_ad                      , STAT=ierr(3) )
271         !
272         DEALLOCATE( rhd_ad      ,                               &
273            &        rhop_ad     ,                               &
274            &        sshb_ad     , sshn_ad   , ssha_ad   ,       &
275            &        sshu_b_ad   , sshu_n_ad , sshu_a_ad ,       &
276            &        sshv_b_ad   , sshv_n_ad , sshv_a_ad ,       &
277            &                      sshf_n_ad ,                   &
278            &        spgu_ad     , spgv_ad   ,                   &
279            &        gtsu_ad     , gtsv_ad,                      &
280            &        gru_ad      , grv_ad                      , STAT=ierr(4) )
281
282#if defined key_zdfddm
283         DEALLOCATE( rrau_ad, STAT=ierr(5) )
284#endif
285
286         ll_allocad = .FALSE.
287      END IF
288      oce_dealloc_tam = MAXVAL( ierr )
289      IF( oce_dealloc_tam /= 0 )   CALL ctl_warn('oce_dealloc_tam: failed to deallocate arrays')
290      !
291   END FUNCTION oce_dealloc_tam
292   !
293   SUBROUTINE oce_tam_init( kmode )
294      !!-----------------------------------------------------------------------
295      !!
296      !!                  ***  ROUTINE oce_tam_init  ***
297      !!
298      !! ** Purpose : Allocate and initialize the tangent linear and
299      !!              adjoint arrays
300      !!
301      !! ** Method  : kindic = 0  allocate/initialize both tl and ad variables
302      !!              kindic = 1  allocate/initialize only tl variables
303      !!              kindic = 2  allocate/initialize only ad variables
304      !!
305      !! ** Action  :
306      !!
307      !! References :
308      !!
309      !! History :
310      !!        ! 2009-03 (A. Weaver) Initial version (based on oce_tam_init)
311      !!        ! 2010-04 (A. Vidard) Nemo3.2 update
312      !!        ! 2012-09 (A. Vidard) Nemo3.4 update
313      !!-----------------------------------------------------------------------
314      !! * Arguments
315      INTEGER, INTENT(IN) :: kmode
316      INTEGER :: ierr
317      IF ( ( kmode == 0 ) .OR. ( kmode == 1 ) ) THEN
318         IF ( .NOT. ll_alloctl ) ierr = oce_alloc_tam ( 1 )
319         ub_tl(:,:,:) = 0._wp
320         vb_tl(:,:,:) = 0._wp
321         un_tl(:,:,:) = 0._wp
322         vn_tl(:,:,:) = 0._wp
323         ua_tl(:,:,:) = 0._wp
324         va_tl(:,:,:) = 0._wp
325         wn_tl(:,:,:) = 0._wp
326         rotb_tl(:,:,:) = 0._wp
327         rotn_tl(:,:,:) = 0._wp
328         hdivb_tl(:,:,:) = 0._wp
329         hdivn_tl(:,:,:) = 0._wp
330         tsb_tl(:,:,:,:) = 0._wp
331         tsn_tl(:,:,:,:) = 0._wp
332         tsa_tl(:,:,:,:) = 0._wp
333         rn2_tl(:,:,:) = 0._wp
334         rhd_tl(:,:,:) = 0._wp
335         rn2b_tl(:,:,:) = 0._wp
336         rhop_tl(:,:,:) = 0._wp
337         sshb_tl(:,:) = 0._wp
338         sshn_tl(:,:) = 0._wp
339         ssha_tl(:,:) = 0._wp
340         sshu_b_tl(:,:) = 0._wp
341         sshu_n_tl(:,:) = 0._wp
342         sshu_a_tl(:,:) = 0._wp
343         sshv_b_tl(:,:) = 0._wp
344         sshv_n_tl(:,:) = 0._wp
345         sshv_a_tl(:,:) = 0._wp
346         sshf_n_tl(:,:) = 0._wp
347         spgu_tl(:,:) = 0._wp
348         spgv_tl(:,:) = 0._wp
349         gtsu_tl(:,:,:) = 0._wp
350         gtsv_tl(:,:,:) = 0._wp
351         gru_tl(:,:) = 0._wp
352         grv_tl(:,:) = 0._wp
353#if defined key_zdfddm
354         rrau_tl(:,:,:) = 0._wp
355#endif
356      END IF
357      IF ( ( kmode == 0 ) .OR. ( kmode == 2 ) ) THEN
358         IF ( .NOT. ll_allocad ) ierr = oce_alloc_tam ( 2 )
359         ub_ad(:,:,:) = 0._wp
360         vb_ad(:,:,:) = 0._wp
361         un_ad(:,:,:) = 0._wp
362         vn_ad(:,:,:) = 0._wp
363         ua_ad(:,:,:) = 0._wp
364         va_ad(:,:,:) = 0._wp
365         wn_ad(:,:,:) = 0._wp
366         rotb_ad(:,:,:) = 0._wp
367         rotn_ad(:,:,:) = 0._wp
368         hdivb_ad(:,:,:) = 0._wp
369         hdivn_ad(:,:,:) = 0._wp
370         tsb_ad(:,:,:,:) = 0._wp
371         tsn_ad(:,:,:,:) = 0._wp
372         tsa_ad(:,:,:,:) = 0._wp
373         rn2_ad(:,:,:) = 0._wp
374         rhd_ad(:,:,:) = 0._wp
375         rn2b_ad(:,:,:) = 0._wp
376         rhop_ad(:,:,:) = 0._wp
377         sshb_ad(:,:) = 0._wp
378         sshn_ad(:,:) = 0._wp
379         ssha_ad(:,:) = 0._wp
380         sshu_b_ad(:,:) = 0._wp
381         sshu_n_ad(:,:) = 0._wp
382         sshu_a_ad(:,:) = 0._wp
383         sshv_b_ad(:,:) = 0._wp
384         sshv_n_ad(:,:) = 0._wp
385         sshv_a_ad(:,:) = 0._wp
386         sshf_n_ad(:,:) = 0._wp
387         spgu_ad(:,:) = 0._wp
388         spgv_ad(:,:) = 0._wp
389         gtsu_ad(:,:,:) = 0._wp
390         gtsv_ad(:,:,:) = 0._wp
391         gru_ad(:,:) = 0._wp
392         grv_ad(:,:) = 0._wp
393         !
394#if defined key_zdfddm
395         rrau_ad(:,:,:) = 0._Wp
396#endif
397      END IF
398
399   END SUBROUTINE oce_tam_init
400
401END MODULE oce_tam
Note: See TracBrowser for help on using the repository browser.