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 branches/2012/dev_r3604_LEGI8_TAM/NEMOGCM/NEMO/OPATAM_SRC – NEMO

source: branches/2012/dev_r3604_LEGI8_TAM/NEMOGCM/NEMO/OPATAM_SRC/oce_tam.F90 @ 3611

Last change on this file since 3611 was 3611, checked in by pabouttier, 11 years ago

Add TAM code and ORCA2_TAM configuration - see Ticket #1007

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