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/TAM_V3_0/NEMOTAM/OPATAM_SRC – NEMO

source: branches/TAM_V3_0/NEMOTAM/OPATAM_SRC/oce_tam.F90 @ 1885

Last change on this file since 1885 was 1885, checked in by rblod, 14 years ago

add TAM sources

File size: 25.5 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   !!----------------------------------------------------------------------
20   !!   oce_tam_init : Allocate and initialize the TAM fields
21   !!----------------------------------------------------------------------
22   !! * Modules used   
23   USE par_kind, ONLY : &  ! Precision variables
24      & wp
25   USE par_oce, ONLY : &   ! Ocean space and time domain variables
26      & jpi, &
27      & jpj, &
28      & jpk
29
30   IMPLICIT NONE
31
32   !! * Routine accessibility
33   PRIVATE
34   PUBLIC &
35      & oce_tam_init, & !: Initialize the TAM fields
36                    !:
37      & ub_tl,    & !: Tangent linear of before u-component velocity
38      & un_tl,    & !: Tangent linear of now u-component velocity
39      & ua_tl,    & !: Tangent linear of after u-component velocity
40      & vb_tl,    & !: Tangent linear of before v-component velocity
41      & vn_tl,    & !: Tangent linear of now v-component velocity
42      & va_tl,    & !: Tangent linear of after v-component velocity
43      & wn_tl,    & !: Tangent linear of now w-component velocity
44      & rotb_tl,  & !: Tangent linear of now relative vorticity
45      & rotn_tl,  & !: Tangent linear of now relative vorticity
46      & hdivb_tl, & !: Tangent linear of before horizontal divergence
47      & hdivn_tl, & !: Tangent linear of now horizontal divergence
48      & tb_tl,    & !: Tangent linear of before of temperature
49      & tn_tl,    & !: Tangent linear of now temperature
50      & ta_tl,    & !: Tangent linear of after temperature
51      & sb_tl,    & !: Tangent linear of before salinity
52      & sn_tl,    & !: Tangent linear of now salinity
53      & sa_tl,    & !: Tangent linear of after salinity
54      & rhd_tl,   & !: Tangent linear of Now in situ density anomaly 
55      & rhop_tl,  & !: Tangent linear of potential volumic mass (kg/m3)
56      & rn2_tl,   & !: Tangent linear of Now Brunt-Vaisala frequency
57                    !:
58      & spgu_tl,  & !: Tangent linear of horizontal surface pressure gradient
59      & spgv_tl,  & !:
60#if defined key_dynspg_rl
61      & bsfb_tl,  & !: Tangent linear of before barotropic streamfunction
62      & bsfn_tl,  & !: Tangent linear of now barotropic streamfunction
63      & bsfd_tl,  & !: Tangent linear of now trend of barotropic streamfunction
64#else
65      & sshb_tl,  & !: Tangent linear of before sea surface height
66      & sshn_tl,  & !: Tangent linear of now sea surface height
67      & ssha_tl,  & !: Tangent linear of after sea surface height
68      & sshu_tl,  & !: Tangent linear of sea surface height at u-point
69      & sshv_tl,  & !: Tangent linear of sea surface height at v-point
70      & sshbb_tl, & !: Tangent linear of before before sea surface height
71#endif
72      & gtu_tl,   & !: Tangent Linear of
73      & gsu_tl,   & !:
74      & gru_tl,   & !: t-, s- and rd horizontal gradient at u- and
75      & gtv_tl,   & !:   v-points at bottom ocean level
76      & gsv_tl,   & !:
77      & grv_tl,   & !:
78                    !:
79      & ub_ad,    & !: Adjoint of before u-component velocity
80      & un_ad,    & !: Adjoint of now u-component velocity
81      & ua_ad,    & !: Adjoint of after u-component velocity
82      & vb_ad,    & !: Adjoint of before v-component velocity
83      & vn_ad,    & !: Adjoint of now v-component velocity
84      & va_ad,    & !: Adjoint of after v-component velocity
85      & wn_ad,    & !: Adjoint of now w-component velocity
86      & rotb_ad,  & !: Adjoint of before relative vorticity
87      & rotn_ad,  & !: Adjoint of now relative vorticity
88      & hdivb_ad, & !: Adjoint of before horizontal divergence
89      & hdivn_ad, & !: Adjoint of now horizontal divergence
90      & tb_ad,    & !: Adjoint of before temperature
91      & tn_ad,    & !: Adjoint of now temperature
92      & ta_ad,    & !: Adjoint of after temperature
93      & sb_ad,    & !: Adjoint of before salinity
94      & sn_ad,    & !: Adjoint of now salinity
95      & sa_ad,    & !: Adjoint of after salinity
96      & rhd_ad,   & !: Adjoint of now in situ density anomaly 
97      & rhop_ad,  & !: Adjoint of potential volumic mass (kg/m3)
98      & rn2_ad,   & !: Adjoint of now Brunt-Vaisala frequency
99      & spgu_ad,  & !: Adjoint linear of horizontal surface pressure gradient
100      & spgv_ad,  & !:
101#if defined key_dynspg_rl
102      & bsfb_ad,  & !: Adjoint of before barotropic streamfunction
103      & bsfn_ad,  & !: Adjoint of now barotropic streamfunction
104      & bsfd_ad,  & !: Adjoint of now trend of barotropic streamfunction
105#else
106      & sshb_ad,  & !: Adjoint of before sea surface height
107      & sshn_ad,  & !: Adjoint of now sea surface height
108      & ssha_ad,  & !: Adjoint of after sea surface height
109      & sshu_ad,  & !: Adjoint of sea surface height at u-point
110      & sshv_ad,  & !: Adjoint of sea surface height at v-point
111      & sshbb_ad, & !: Adjoint of before before sea surface height
112#endif
113      & gtu_ad,   & !: Adjoint of
114      & gsu_ad,   & !:
115      & gru_ad,   & !: t-, s- and rd horizontal gradient at u- and
116      & gtv_ad,   & !:   v-points at bottom ocean level
117      & gsv_ad,   & !:
118      & grv_ad
119
120#if defined key_zdfddm
121!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
122!!!!  AW: The declaration/allocation/initialization of these variables
123!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
124!!!!      with NEMO.
125   PUBLIC &
126      & rrau_tl,  & !: Tangent Linear of heat/salt buoyancy flux ratio
127      & rrau_ad     !: Adjoint of heat/salt buoyancy flux ratio
128!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
129#endif
130
131   !! * Module variables
132   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
133      & ub_tl,    & !: Tangent linear of before u-component velocity
134      & un_tl,    & !: Tangent linear of now u-component velocity
135      & ua_tl,    & !: Tangent linear of after u-component velocity
136      & vb_tl,    & !: Tangent linear of before v-component velocity
137      & vn_tl,    & !: Tangent linear of now v-component velocity
138      & va_tl,    & !: Tangent linear of after v-component velocity
139      & wn_tl,    & !: Tangent linear of now w-component velocity
140      & rotb_tl,  & !: Tangent linear of before relative vorticity
141      & rotn_tl,  & !: Tangent linear of now relative vorticity
142      & hdivb_tl, & !: Tangent linear of before horizontal divergence
143      & hdivn_tl, & !: Tangent linear of now horizontal divergence
144      & tb_tl,    & !: Tangent linear of before temperature
145      & tn_tl,    & !: Tangent linear of now temperature
146      & ta_tl,    & !: Tangent linear of after temperature
147      & sb_tl,    & !: Tangent linear of before salinity
148      & sn_tl,    & !: Tangent linear of now salinity
149      & sa_tl,    & !: Tangent linear of after salinity
150      & rhd_tl,   & !: Tangent linear of now in situ density anomaly 
151      & rhop_tl,  & !: Tangent linear of potential volumic mass (kg/m3)
152      & rn2_tl      !: Tangent linear of now Brunt-Vaisala frequency
153
154   REAL(KIND=wp), DIMENSION(:,:), ALLOCATABLE :: &
155      & spgu_tl,  & !: Tangent of horizontal surface pressure gradient
156      & spgv_tl,  & !:
157#if defined key_dynspg_rl
158      & bsfb_tl,  & !: Tangent linear of before barotropic streamfunction
159      & bsfn_tl,  & !: Tangent linear of now barotropic streamfunction
160      & bsfd_tl,  & !: Tangent linear of now trend of barotropic streamfunction
161#else 
162      & sshb_tl,  & !: Tangent linear of before sea surface height
163      & sshn_tl,  & !: Tangent linear of now sea surface height
164      & ssha_tl,  & !: Tangent linear of after sea surface height
165      & sshu_tl,  & !: Tangent linear of sea surface height at u-point
166      & sshv_tl,  & !: Tangent linear of sea surface height at v-point
167      & sshbb_tl, & !: Tangent linear of before before sea surface height
168#endif
169      & gtu_tl,   & !: Tangent Linear of
170      & gsu_tl,   & !:
171      & gru_tl,   & !: t-, s- and rd horizontal gradient at u- and
172      & gtv_tl,   & !:   v-points at bottom ocean level
173      & gsv_tl,   & !:
174      & grv_tl       
175
176   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
177      & ub_ad,    & !: Adjoint of before u-component velocity
178      & un_ad,    & !: Adjoint of now u-component velocity
179      & ua_ad,    & !: Adjoint of after u-component velocity
180      & vb_ad,    & !: Adjoint of before v-component velocity
181      & vn_ad,    & !: Adjoint of now v-component velocity
182      & va_ad,    & !: Adjoint of after v-component velocity
183      & wn_ad,    & !: Adjoint of now w-component velocity
184      & rotb_ad,  & !: Adjoint of before relative vorticity
185      & rotn_ad,  & !: Adjoint of now relative vorticity
186      & hdivb_ad, & !: Adjoint of before horizontal divergence
187      & hdivn_ad, & !: Adjoint of now horizontal divergence
188      & tb_ad,    & !: Adjoint of before temperature
189      & tn_ad,    & !: Adjoint of now temperature
190      & ta_ad,    & !: Adjoint of after temperature
191      & sb_ad,    & !: Adjoint of before salinity
192      & sn_ad,    & !: Adjoint of now salinity
193      & sa_ad,    & !: Adjoint of after salinity
194      & rhd_ad,   & !: Adjoint of now in situ density anomaly 
195      & rhop_ad,  & !: Adjoint of potential volumic mass (kg/m3)
196      & rn2_ad      !: Adjoint of now Brunt-Vaisala frequency
197
198#if defined key_zdfddm
199!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
200!!!!  AW: The declaration/allocation/initialization of these variables
201!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
202!!!!      with NEMO.
203   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
204      & rrau_tl,  & !: Tangent Linear of heat/salt buoyancy flux ratio
205      & rrau_ad     !: Adjoint of heat/salt buoyancy flux ratio
206!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
207#endif
208
209   REAL(KIND=wp), DIMENSION(:,:), ALLOCATABLE :: &
210      & spgu_ad,  & !: Adjoint of horizontal surface pressure gradient
211      & spgv_ad,  & !:
212#if defined key_dynspg_rl
213      & bsfb_ad,  & !: Adjoint of before barotropic streamfunction
214      & bsfn_ad,  & !: Adjoint of now barotropic streamfunction
215      & bsfd_ad,  & !: Adjoint of now trend of barotropic streamfunction
216#else 
217      & sshb_ad,  & !: Adjoint of before sea surface height
218      & sshn_ad,  & !: Adjoint of now sea surface height
219      & ssha_ad,  & !: Adjoint of after sea surface height
220      & sshu_ad,  & !: Adjoint of sea surface height at u-point
221      & sshv_ad,  & !: Adjoint of sea surface height at v-point
222      & sshbb_ad, & !: Adjoint of before before sea surface height
223#endif
224      & gtu_ad,   & !: Adjoint of
225      & gsu_ad,   & !:
226      & gru_ad,   & !: t-, s- and rd horizontal gradient at u- and
227      & gtv_ad,   & !:   v-points at bottom ocean level
228      & gsv_ad,   & !:
229      & grv_ad
230
231CONTAINS
232
233   SUBROUTINE oce_tam_init( kindic )
234      !!-----------------------------------------------------------------------
235      !!
236      !!                  ***  ROUTINE oce_tam_init  ***
237      !!
238      !! ** Purpose : Allocate and initialize the tangent linear and
239      !!              adjoint arrays
240      !!
241      !! ** Method  : kindic = 0  allocate/initialize both tl and ad variables
242      !!              kindic = 1  allocate/initialize only tl variables
243      !!              kindic = 2  allocate/initialize only ad variables
244      !!
245      !! ** Action  :
246      !!                   
247      !! References :
248      !!
249      !! History :
250      !!        ! 07-07 (K. Mogensen) Initial version
251      !!        ! 08-07 (A. Vidard) Add missing variables
252      !!        ! 08-07 (A. Weaver) Add missing variables
253      !!        ! 09-03 (A. Weaver) Combine tl_init/ad_init into oce_tam_init,
254      !!                            make consistent with oce.F90 in NEMO
255      !!        ! 09-03 (A. Vidard) bug fixes on ssh*_ad alloc
256      !!-----------------------------------------------------------------------
257      !! * Arguments
258      INTEGER, INTENT(IN) :: &
259         & kindic        ! indicate which variables to allocate/initialize
260
261      !! * Local declarations
262      ! Allocate tangent linear variable arrays
263      ! ---------------------------------------
264     
265      IF ( kindic == 0 .OR. kindic == 1 ) THEN
266
267         IF ( .NOT. ALLOCATED(ub_tl) ) THEN
268
269            ALLOCATE( ub_tl(jpi,jpj,jpk) )
270
271         ENDIF
272         IF ( .NOT. ALLOCATED(un_tl) ) THEN
273
274            ALLOCATE( un_tl(jpi,jpj,jpk) )
275
276         ENDIF
277         IF ( .NOT. ALLOCATED(ua_tl) ) THEN
278
279            ALLOCATE( ua_tl(jpi,jpj,jpk) )
280
281         ENDIF
282         IF ( .NOT. ALLOCATED(vb_tl) ) THEN
283
284            ALLOCATE( vb_tl(jpi,jpj,jpk) )
285
286         ENDIF
287         IF ( .NOT. ALLOCATED(vn_tl) ) THEN
288
289            ALLOCATE( vn_tl(jpi,jpj,jpk) )
290
291         ENDIF
292         IF ( .NOT. ALLOCATED(va_tl) ) THEN
293
294            ALLOCATE( va_tl(jpi,jpj,jpk) )
295
296         ENDIF
297         IF ( .NOT. ALLOCATED(wn_tl) ) THEN
298
299            ALLOCATE( wn_tl(jpi,jpj,jpk) )
300           
301         ENDIF
302         IF ( .NOT. ALLOCATED(rotb_tl) ) THEN
303
304            ALLOCATE( rotb_tl(jpi,jpj,jpk) )
305           
306         ENDIF
307         IF ( .NOT. ALLOCATED(rotn_tl) ) THEN
308
309            ALLOCATE( rotn_tl(jpi,jpj,jpk) )
310           
311         ENDIF
312         IF ( .NOT. ALLOCATED(hdivb_tl) ) THEN
313
314            ALLOCATE( hdivb_tl(jpi,jpj,jpk) )
315           
316         ENDIF
317         IF ( .NOT. ALLOCATED(hdivn_tl) ) THEN
318
319            ALLOCATE( hdivn_tl(jpi,jpj,jpk) )
320           
321         ENDIF
322         IF ( .NOT. ALLOCATED(tb_tl) ) THEN
323
324            ALLOCATE( tb_tl(jpi,jpj,jpk) )
325           
326         ENDIF
327         IF ( .NOT. ALLOCATED(tn_tl) ) THEN
328
329            ALLOCATE( tn_tl(jpi,jpj,jpk) )
330           
331         ENDIF
332         IF ( .NOT. ALLOCATED(ta_tl) ) THEN
333
334            ALLOCATE( ta_tl(jpi,jpj,jpk) )
335           
336         ENDIF
337         IF ( .NOT. ALLOCATED(sb_tl) ) THEN
338
339            ALLOCATE( sb_tl(jpi,jpj,jpk) )
340           
341         ENDIF
342         IF ( .NOT. ALLOCATED(sn_tl) ) THEN
343
344            ALLOCATE( sn_tl(jpi,jpj,jpk) )
345           
346         ENDIF
347         IF ( .NOT. ALLOCATED(sa_tl) ) THEN
348
349            ALLOCATE( sa_tl(jpi,jpj,jpk) )
350           
351         ENDIF
352         IF ( .NOT. ALLOCATED(rhd_tl) ) THEN
353
354            ALLOCATE( rhd_tl(jpi,jpj,jpk) )
355           
356         ENDIF
357         IF ( .NOT. ALLOCATED(rhop_tl) ) THEN
358
359            ALLOCATE( rhop_tl(jpi,jpj,jpk) )
360           
361         ENDIF
362         IF ( .NOT. ALLOCATED(rn2_tl) ) THEN
363
364            ALLOCATE( rn2_tl(jpi,jpj,jpk) )
365
366         ENDIF
367         IF ( .NOT. ALLOCATED(spgu_tl) ) THEN
368
369            ALLOCATE( spgu_tl(jpi,jpj) )
370           
371         ENDIF
372         IF ( .NOT. ALLOCATED(spgv_tl) ) THEN
373
374            ALLOCATE( spgv_tl(jpi,jpj) )
375           
376         ENDIF
377#if defined key_dynspg_rl
378         IF ( .NOT. ALLOCATED(bsfb_tl) ) THEN
379
380            ALLOCATE( bsfb_tl(jpi,jpj) )
381         
382         ENDIF
383         IF ( .NOT. ALLOCATED(bsfn_tl) ) THEN
384
385            ALLOCATE( bsfn_tl(jpi,jpj) )
386         
387         ENDIF
388         IF ( .NOT. ALLOCATED(bsfd_tl) ) THEN
389
390            ALLOCATE( bsfd_tl(jpi,jpj) )
391         
392         ENDIF
393#else
394         IF (.NOT. ALLOCATED(sshb_tl) ) THEN
395         
396            ALLOCATE( sshb_tl(jpi,jpj) )
397
398         ENDIF
399         IF (.NOT. ALLOCATED(sshn_tl) ) THEN
400         
401            ALLOCATE( sshn_tl(jpi,jpj) )
402
403         ENDIF
404         IF (.NOT. ALLOCATED(ssha_tl) ) THEN
405         
406            ALLOCATE( ssha_tl(jpi,jpj) )
407
408         ENDIF
409         IF (.NOT. ALLOCATED(sshu_tl) ) THEN
410         
411            ALLOCATE( sshu_tl(jpi,jpj) )
412
413         ENDIF
414         IF (.NOT. ALLOCATED(sshv_tl) ) THEN
415         
416            ALLOCATE( sshv_tl(jpi,jpj) )
417
418         ENDIF
419         IF (.NOT. ALLOCATED(sshbb_tl) ) THEN
420         
421            ALLOCATE( sshbb_tl(jpi,jpj) )
422
423         ENDIF
424#endif     
425         IF ( .NOT. ALLOCATED(gtu_tl) ) THEN
426
427            ALLOCATE( gtu_tl(jpi,jpj) )
428
429         ENDIF
430         IF ( .NOT. ALLOCATED(gtv_tl) ) THEN
431
432            ALLOCATE( gtv_tl(jpi,jpj) )
433
434         ENDIF
435         IF ( .NOT. ALLOCATED(gsu_tl) ) THEN
436
437            ALLOCATE( gsu_tl(jpi,jpj) )
438
439         ENDIF
440         IF ( .NOT. ALLOCATED(gsv_tl) ) THEN
441
442            ALLOCATE( gsv_tl(jpi,jpj) )
443
444         ENDIF
445         IF ( .NOT. ALLOCATED(gru_tl) ) THEN
446
447            ALLOCATE( gru_tl(jpi,jpj) )
448
449         ENDIF
450         IF ( .NOT. ALLOCATED(grv_tl) ) THEN
451
452            ALLOCATE( grv_tl(jpi,jpj) )
453
454         ENDIF
455
456#if defined key_zdfddm
457!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
458!!!!  AW: The declaration/allocation/initialization of these variables
459!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
460!!!!      with NEMO.
461         IF ( .NOT. ALLOCATED(rrau_tl) ) THEN
462
463            ALLOCATE( rrau_tl(jpi,jpj,jpk) )
464           
465         ENDIF
466!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
467#endif
468
469         ! Initialize tangent linear variable arrays to zero
470         ! -------------------------------------------------
471
472         ub_tl   (:,:,:) = 0.0_wp
473         un_tl   (:,:,:) = 0.0_wp
474         ua_tl   (:,:,:) = 0.0_wp
475         vb_tl   (:,:,:) = 0.0_wp
476         vn_tl   (:,:,:) = 0.0_wp
477         va_tl   (:,:,:) = 0.0_wp
478         wn_tl   (:,:,:) = 0.0_wp
479         rotb_tl (:,:,:) = 0.0_wp
480         rotn_tl (:,:,:) = 0.0_wp
481         hdivb_tl(:,:,:) = 0.0_wp
482         hdivn_tl(:,:,:) = 0.0_wp
483         tb_tl   (:,:,:) = 0.0_wp
484         tn_tl   (:,:,:) = 0.0_wp
485         ta_tl   (:,:,:) = 0.0_wp
486         sb_tl   (:,:,:) = 0.0_wp
487         sn_tl   (:,:,:) = 0.0_wp
488         sa_tl   (:,:,:) = 0.0_wp
489         rhd_tl  (:,:,:) = 0.0_wp
490         rhop_tl (:,:,:) = 0.0_wp
491         rn2_tl  (:,:,:) = 0.0_wp
492
493         spgu_tl(:,:)  = 0.0_wp
494         spgv_tl(:,:)  = 0.0_wp
495
496#if defined key_dynspg_rl
497         bsfb_tl(:,:)  = 0.0_wp
498         bsfn_tl(:,:)  = 0.0_wp
499         bsfd_tl(:,:)  = 0.0_wp
500#else
501         sshb_tl(:,:)  = 0.0_wp
502         sshn_tl(:,:)  = 0.0_wp
503         ssha_tl(:,:)  = 0.0_wp
504         sshu_tl(:,:)  = 0.0_wp
505         sshv_tl(:,:)  = 0.0_wp
506         sshbb_tl(:,:) = 0.0_wp
507#endif
508         gtu_tl (:,:)  = 0.0_wp
509         gsu_tl (:,:)  = 0.0_wp
510         gru_tl (:,:)  = 0.0_wp
511         gtv_tl (:,:)  = 0.0_wp
512         gsv_tl (:,:)  = 0.0_wp
513         grv_tl (:,:)  = 0.0_wp
514
515#if defined key_zdfddm
516!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
517!!!!  AW: The declaration/allocation/initialization of these variables
518!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
519!!!!      with NEMO.
520         rrau_tl(:,:,:) = 0.0_wp
521!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
522#endif
523
524      ENDIF
525
526      IF ( kindic == 0 .OR. kindic == 2 ) THEN
527
528         ! Allocate adjoint variable arrays
529         ! --------------------------------
530     
531         IF ( .NOT. ALLOCATED(ub_ad) ) THEN
532
533            ALLOCATE( ub_ad(jpi,jpj,jpk) )
534
535         ENDIF
536         IF ( .NOT. ALLOCATED(un_ad) ) THEN
537
538            ALLOCATE( un_ad(jpi,jpj,jpk) )
539
540         ENDIF
541         IF ( .NOT. ALLOCATED(ua_ad) ) THEN
542
543            ALLOCATE( ua_ad(jpi,jpj,jpk) )
544
545         ENDIF
546         IF ( .NOT. ALLOCATED(vb_ad) ) THEN
547
548            ALLOCATE( vb_ad(jpi,jpj,jpk) )
549
550         ENDIF
551         IF ( .NOT. ALLOCATED(vn_ad) ) THEN
552
553            ALLOCATE( vn_ad(jpi,jpj,jpk) )
554
555         ENDIF
556         IF ( .NOT. ALLOCATED(va_ad) ) THEN
557
558            ALLOCATE( va_ad(jpi,jpj,jpk) )
559
560         ENDIF
561         IF ( .NOT. ALLOCATED(wn_ad) ) THEN
562
563            ALLOCATE( wn_ad(jpi,jpj,jpk) )
564           
565         ENDIF
566         IF ( .NOT. ALLOCATED(rotb_ad) ) THEN
567
568            ALLOCATE( rotb_ad(jpi,jpj,jpk) )
569           
570         ENDIF
571         IF ( .NOT. ALLOCATED(rotn_ad) ) THEN
572
573            ALLOCATE( rotn_ad(jpi,jpj,jpk) )
574           
575         ENDIF
576         IF ( .NOT. ALLOCATED(hdivb_ad) ) THEN
577
578            ALLOCATE( hdivb_ad(jpi,jpj,jpk) )
579           
580         ENDIF
581         IF ( .NOT. ALLOCATED(hdivn_ad) ) THEN
582
583            ALLOCATE( hdivn_ad(jpi,jpj,jpk) )
584           
585         ENDIF
586         IF ( .NOT. ALLOCATED(tb_ad) ) THEN
587
588            ALLOCATE( tb_ad(jpi,jpj,jpk) )
589           
590         ENDIF
591         IF ( .NOT. ALLOCATED(tn_ad) ) THEN
592
593            ALLOCATE( tn_ad(jpi,jpj,jpk) )
594           
595         ENDIF
596         IF ( .NOT. ALLOCATED(ta_ad) ) THEN
597
598            ALLOCATE( ta_ad(jpi,jpj,jpk) )
599           
600         ENDIF
601         IF ( .NOT. ALLOCATED(sb_ad) ) THEN
602
603            ALLOCATE( sb_ad(jpi,jpj,jpk) )
604           
605         ENDIF
606         IF ( .NOT. ALLOCATED(sn_ad) ) THEN
607
608            ALLOCATE( sn_ad(jpi,jpj,jpk) )
609           
610         ENDIF
611         IF ( .NOT. ALLOCATED(sa_ad) ) THEN
612
613            ALLOCATE( sa_ad(jpi,jpj,jpk) )
614           
615         ENDIF
616         IF ( .NOT. ALLOCATED(rhd_ad) ) THEN
617
618            ALLOCATE( rhd_ad(jpi,jpj,jpk) )
619           
620         ENDIF
621         IF ( .NOT. ALLOCATED(rhop_ad) ) THEN
622
623            ALLOCATE( rhop_ad(jpi,jpj,jpk) )
624           
625         ENDIF
626         IF ( .NOT. ALLOCATED(rn2_ad) ) THEN
627
628            ALLOCATE( rn2_ad(jpi,jpj,jpk) )
629
630         ENDIF
631         IF ( .NOT. ALLOCATED(spgu_ad) ) THEN
632
633            ALLOCATE( spgu_ad(jpi,jpj) )
634           
635         ENDIF
636         IF ( .NOT. ALLOCATED(spgv_ad) ) THEN
637
638            ALLOCATE( spgv_ad(jpi,jpj) )
639           
640         ENDIF
641#if defined key_dynspg_rl
642         IF ( .NOT. ALLOCATED(bsfb_ad) ) THEN
643
644            ALLOCATE( bsfb_ad(jpi,jpj) )
645         
646         ENDIF
647         IF ( .NOT. ALLOCATED(bsfn_ad) ) THEN
648
649            ALLOCATE( bsfn_ad(jpi,jpj) )
650         
651         ENDIF
652         IF ( .NOT. ALLOCATED(bsfd_ad) ) THEN
653
654            ALLOCATE( bsfd_ad(jpi,jpj) )
655         
656         ENDIF
657#else
658         IF (.NOT. ALLOCATED(sshb_ad) ) THEN
659         
660            ALLOCATE( sshb_ad(jpi,jpj) )
661
662         ENDIF
663         IF (.NOT. ALLOCATED(sshn_ad) ) THEN
664         
665            ALLOCATE( sshn_ad(jpi,jpj) )
666
667         ENDIF
668         IF (.NOT. ALLOCATED(ssha_ad) ) THEN
669         
670            ALLOCATE( ssha_ad(jpi,jpj) )
671
672         ENDIF
673         IF (.NOT. ALLOCATED(sshu_ad) ) THEN
674         
675            ALLOCATE( sshu_ad(jpi,jpj) )
676
677         ENDIF
678         IF (.NOT. ALLOCATED(sshv_ad) ) THEN
679         
680            ALLOCATE( sshv_ad(jpi,jpj) )
681
682         ENDIF
683         IF (.NOT. ALLOCATED(sshbb_ad) ) THEN
684         
685            ALLOCATE( sshbb_ad(jpi,jpj) )
686
687         ENDIF
688#endif     
689         IF ( .NOT. ALLOCATED(gtu_ad) ) THEN
690
691            ALLOCATE( gtu_ad(jpi,jpj) )
692
693         ENDIF
694         IF ( .NOT. ALLOCATED(gtv_ad) ) THEN
695
696            ALLOCATE( gtv_ad(jpi,jpj) )
697
698         ENDIF
699         IF ( .NOT. ALLOCATED(gsu_ad) ) THEN
700
701            ALLOCATE( gsu_ad(jpi,jpj) )
702
703         ENDIF
704         IF ( .NOT. ALLOCATED(gsv_ad) ) THEN
705
706            ALLOCATE( gsv_ad(jpi,jpj) )
707
708         ENDIF
709         IF ( .NOT. ALLOCATED(gru_ad) ) THEN
710
711            ALLOCATE( gru_ad(jpi,jpj) )
712
713         ENDIF
714         IF ( .NOT. ALLOCATED(grv_ad) ) THEN
715
716            ALLOCATE( grv_ad(jpi,jpj) )
717
718         ENDIF
719
720#if defined key_zdfddm
721!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
722!!!!  AW: The declaration/allocation/initialization of these variables
723!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
724!!!!      with NEMO.
725         IF ( .NOT. ALLOCATED(rrau_ad) ) THEN
726
727            ALLOCATE( rrau_ad(jpi,jpj,jpk) )
728           
729         ENDIF
730!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
731#endif
732
733         ! Initialize adjoint variable arrays to zero
734         ! ------------------------------------------
735
736         ub_ad   (:,:,:) = 0.0_wp
737         un_ad   (:,:,:) = 0.0_wp
738         ua_ad   (:,:,:) = 0.0_wp
739         vb_ad   (:,:,:) = 0.0_wp
740         vn_ad   (:,:,:) = 0.0_wp
741         va_ad   (:,:,:) = 0.0_wp
742         wn_ad   (:,:,:) = 0.0_wp
743         rotb_ad (:,:,:) = 0.0_wp
744         rotn_ad (:,:,:) = 0.0_wp
745         hdivb_ad(:,:,:) = 0.0_wp
746         hdivn_ad(:,:,:) = 0.0_wp
747         tb_ad   (:,:,:) = 0.0_wp
748         tn_ad   (:,:,:) = 0.0_wp
749         ta_ad   (:,:,:) = 0.0_wp
750         sb_ad   (:,:,:) = 0.0_wp
751         sn_ad   (:,:,:) = 0.0_wp
752         sa_ad   (:,:,:) = 0.0_wp
753         rhd_ad  (:,:,:) = 0.0_wp
754         rhop_ad (:,:,:) = 0.0_wp
755         rn2_ad  (:,:,:) = 0.0_wp
756
757         spgu_ad(:,:) = 0.0_wp
758         spgv_ad(:,:) = 0.0_wp
759#if defined key_dynspg_rl
760         bsfn_ad(:,:) = 0.0_wp
761         bsfb_ad(:,:) = 0.0_wp
762         bsfd_ad(:,:) = 0.0_wp
763#else
764         sshb_ad(:,:)  = 0.0_wp
765         sshn_ad(:,:)  = 0.0_wp
766         ssha_ad(:,:)  = 0.0_wp
767         sshu_ad(:,:)  = 0.0_wp
768         sshv_ad(:,:)  = 0.0_wp
769         sshbb_ad(:,:) = 0.0_wp
770#endif
771         gtu_ad (:,:) = 0.0_wp
772         gsu_ad (:,:) = 0.0_wp
773         gru_ad (:,:) = 0.0_wp
774         gtv_ad (:,:) = 0.0_wp
775         gsv_ad (:,:) = 0.0_wp
776         grv_ad (:,:) = 0.0_wp
777
778#if defined key_zdfddm
779!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
780!!!!  AW: The declaration/allocation/initialization of these variables
781!!!!      should be moved to a new module zdf_ddm_tam_init to be consistent
782!!!!      with NEMO.
783         rrau_ad  (:,:,:) = 0.0_wp
784!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
785#endif
786
787      ENDIF
788
789   END SUBROUTINE oce_tam_init
790     
791END MODULE oce_tam
Note: See TracBrowser for help on using the repository browser.