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.
tranxt_tam.F90 in branches/TAM_V3_0/NEMOTAM/OPATAM_SRC/TRA – NEMO

source: branches/TAM_V3_0/NEMOTAM/OPATAM_SRC/TRA/tranxt_tam.F90 @ 2587

Last change on this file since 2587 was 2587, checked in by vidard, 13 years ago

refer to ticket #798

  • Property svn:executable set to *
File size: 29.8 KB
Line 
1MODULE tranxt_tam
2#ifdef key_tam
3   !!======================================================================
4   !!                       ***  MODULE  tranxt_tam  ***
5   !! Ocean active tracers:  time stepping on temperature and salinity
6   !!               Tangent and Adjoint module
7   !!======================================================================
8   !! History of the direct module:
9   !!            7.0  !  91-11  (G. Madec)  Original code
10   !!                 !  93-03  (M. Guyon)  symetrical conditions
11   !!                 !  96-02  (G. Madec & M. Imbard)  opa release 8.0
12   !!            8.0  !  96-04  (A. Weaver)  Euler forward step
13   !!            8.2  !  99-02  (G. Madec, N. Grima)  semi-implicit pressure grad.
14   !!  NEMO      1.0  !  2002-08  (G. Madec)  F90: Free form and module
15   !!             -   !  2002-11  (C. Talandier, A-M Treguier) Open boundaries
16   !!             -   !  2005-04  (C. Deltel) Add Asselin trend in the ML budget
17   !!            2.0  !  2006-02  (L. Debreu, C. Mazauric) Agrif implementation
18   !!            3.0  !  2008-06  (G. Madec)  time stepping always done in trazd   !! History of the TAM module:
19   !!            2.0  !  2008-09  (A. Vidard) tam of the 2006-02 version
20   !!            3.0  !  2008-11  (A. Vidard) tam of the 2008-06 version
21   !!             -   !  2009-01  (A. Weaver) corrections to test
22   !!----------------------------------------------------------------------
23
24   !!----------------------------------------------------------------------
25   !!   tra_nxt_tan    : time stepping on temperature and salinity (tangent)
26   !!   tra_nxt_adj    : time stepping on temperature and salinity (adjoint)
27   !!----------------------------------------------------------------------
28   USE par_kind      , ONLY: & ! Precision variables
29      & wp
30   USE par_oce       , ONLY: & ! Ocean space and time domain variables
31      & jpi,                 &
32      & jpj,                 & 
33      & jpk,                 &
34      & jpkm1,               &
35      & jpiglo
36   USE oce           , ONLY: &! ocean dynamics and tracers variables
37      & ln_dynhpg_imp
38   USE oce_tam       , ONLY: &! ocean dynamics and tracers variables
39      & tn_tl,               &
40      & tb_tl,               &
41      & ta_tl,               &
42      & sn_tl,               &
43      & sb_tl,               &
44      & sa_tl,               &
45      & tn_ad,               &
46      & tb_ad,               &
47      & ta_ad,               &
48      & sn_ad,               &
49      & sb_ad,               &
50      & sa_ad
51   USE zdf_oce       , ONLY: &
52      & ln_zdfexp
53   USE dom_oce       , ONLY: & ! ocean space and time domain variables
54      & neuler,              &
55      & rdt,                 & 
56      & atfp,                & 
57      & atfp1,               & 
58      & e1t,                 &
59      & e2t,                 &
60# if defined key_vvl
61      & e3t_1,               &
62# else
63#  if defined key_zco
64      & e3t_0,               &
65#  else
66      & e3t,                 &
67#  endif
68# endif
69      & tmask,               &
70      & mig,                 &
71      & mjg,                 &
72      & nldi,                &
73      & nldj,                &
74      & nlei,                &
75      & nlej
76   USE in_out_manager, ONLY: & ! I/O manager
77      & lwp,                 &
78      & numout,              &
79      & nitend,              &
80      & nit000
81   USE lbclnk        , ONLY: &
82      & lbc_lnk
83   USE lbclnk_tam    , ONLY: &
84      & lbc_lnk_adj
85   USE gridrandom    , ONLY: & ! Random Gaussian noise on grids
86      & grid_random
87   USE dotprodfld    , ONLY: & ! Computes dot product for 3D and 2D fields
88      & dot_product
89   USE paresp        , ONLY: & ! Weights for an energy-type scalar product
90      & wesp_t,              &
91      & wesp_s
92   USE tstool_tam    , ONLY: &
93      & prntst_adj,          & !
94      & stdt,                & ! stdev for temperature
95      & stds                   !           salinity
96#if defined key_obc
97#  if defined key_pomme_r025
98   USE obc_oce
99   USE obctra_tam
100#  else
101   Error, OBC not ready.
102#  endif
103#endif
104
105   IMPLICIT NONE
106   PRIVATE
107
108   !! * Routine accessibility
109   PUBLIC   tra_nxt_tan     ! routine called by step_tam.F90
110   PUBLIC   tra_nxt_adj     ! routine called by step_tam.F90
111   PUBLIC   tra_nxt_adj_tst ! routine called by tst.F90
112
113   !! * Substitutions
114#  include "domzgr_substitute.h90"
115
116CONTAINS
117
118   SUBROUTINE tra_nxt_tan( kt )
119      !!----------------------------------------------------------------------
120      !!                   ***  ROUTINE tranxt_tan  ***
121      !!
122      !! ** Purpose of the direct routine:
123      !!             Apply the boundary condition on the after temperature 
124      !!             and salinity fields, achieved the time stepping by adding
125      !!             the Asselin filter on now fields and swapping the fields.
126      !!
127      !! ** Method  :   At this stage of the computation, ta and sa are the
128      !!             after temperature and salinity as the time stepping has
129      !!             been performed in trazdf_imp or trazdf_exp module.
130      !!
131      !!              - Apply lateral boundary conditions on (ta,sa)
132      !!             at the local domain   boundaries through lbc_lnk call,
133      !!             at the radiative open boundaries (lk_obc=T),
134      !!             at the relaxed   open boundaries (lk_bdy=T), and
135      !!             at the AGRIF zoom     boundaries (lk_agrif=T)
136      !!
137      !!              - Apply the Asselin time filter on now fields,
138      !!             save in (ta,sa) an average over the three time levels
139      !!             which will be used to compute rdn and thus the semi-implicit
140      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T), and
141      !!             swap tracer fields to prepare the next time_step.
142      !!                This can be summurized for tempearture as:
143      !!                    zt = (ta+2tn+tb)/4       ln_dynhpg_imp = T
144      !!                    zt = 0                   otherwise
145      !!                    tb = tn + atfp*[ tb - 2 tn + ta ]
146      !!                    tn = ta
147      !!                    ta = zt        (NB: reset to 0 after eos_bn2 call)
148      !!
149      !! ** Action  : - update (tb,sb) and (tn,sn)
150      !!              - (ta,sa) time averaged (t,s)      (ln_dynhpg_imp = T)
151      !!----------------------------------------------------------------------
152      !!
153      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
154      !!
155      !!
156      INTEGER  ::   ji, jj, jk    ! dummy loop indices
157      REAL(wp) ::   zttl, zstl    ! temporary scalars
158      REAL(wp) ::   zfact         ! temporary scalar
159
160      !!----------------------------------------------------------------------
161      IF( kt == nit000 ) THEN
162         IF(lwp) WRITE(numout,*)
163         IF(lwp) WRITE(numout,*) 'tra_nxt_tan : achieve the time stepping by Asselin filter and array swap'
164         IF(lwp) WRITE(numout,*) '~~~~~~~'
165      ENDIF
166
167      ! Update after tracer on domain lateral boundaries
168      !
169      CALL lbc_lnk( ta_tl, 'T', 1. )      ! local domain boundaries  (T-point, unchanged sign)
170      CALL lbc_lnk( sa_tl, 'T', 1. )
171      !
172#if defined key_obc
173      CALL obc_tra_tan( kt )               ! OBC open boundaries
174#endif
175#if defined key_bdy
176      error "bdy not available in tangent yet"
177#endif
178#if defined key_agrif
179      error "agrif not available in tangent yet"
180#endif
181
182      ! Asselin time filter and swap of arrays
183      !                                             
184      IF( neuler == 0 .AND. kt == nit000 ) THEN      ! Euler 1st time step : swap only
185         DO jk = 1, jpkm1
186            tb_tl(:,:,jk) = tn_tl(:,:,jk)                       ! ta, sa remain at their values which
187            sb_tl(:,:,jk) = sn_tl(:,:,jk)                       ! correspond to tn, sn after the sawp
188            tn_tl(:,:,jk) = ta_tl(:,:,jk)                     
189            sn_tl(:,:,jk) = sa_tl(:,:,jk)
190         END DO
191         !                                           
192      ELSE                                           ! Leap-Frog : filter + swap
193         !                                           
194         IF( ln_dynhpg_imp ) THEN                         ! semi-implicit hpg case
195            DO jk = 1, jpkm1                              ! (save the averaged of the 3 time steps
196               DO jj = 1, jpj                             !                   in the after fields)
197                  DO ji = 1, jpi
198                     zttl = ( ta_tl(ji,jj,jk) + 2. * tn_tl(ji,jj,jk) + tb_tl(ji,jj,jk) ) * 0.25
199                     zstl = ( sa_tl(ji,jj,jk) + 2. * sn_tl(ji,jj,jk) + sb_tl(ji,jj,jk) ) * 0.25
200                     tb_tl(ji,jj,jk) = atfp * ( tb_tl(ji,jj,jk) + ta_tl(ji,jj,jk) ) + atfp1 * tn_tl(ji,jj,jk)
201                     sb_tl(ji,jj,jk) = atfp * ( sb_tl(ji,jj,jk) + sa_tl(ji,jj,jk) ) + atfp1 * sn_tl(ji,jj,jk)
202                     tn_tl(ji,jj,jk) = ta_tl(ji,jj,jk)
203                     sn_tl(ji,jj,jk) = sa_tl(ji,jj,jk)
204                     ta_tl(ji,jj,jk) = zttl
205                     sa_tl(ji,jj,jk) = zstl
206                  END DO
207               END DO
208            END DO
209         ELSE                                            ! explicit hpg case
210            DO jk = 1, jpkm1
211               DO jj = 1, jpj
212                  DO ji = 1, jpi
213                     tb_tl(ji,jj,jk) = atfp * ( tb_tl(ji,jj,jk) + ta_tl(ji,jj,jk) ) + atfp1 * tn_tl(ji,jj,jk)
214                     sb_tl(ji,jj,jk) = atfp * ( sb_tl(ji,jj,jk) + sa_tl(ji,jj,jk) ) + atfp1 * sn_tl(ji,jj,jk)
215                     tn_tl(ji,jj,jk) = ta_tl(ji,jj,jk)
216                     sn_tl(ji,jj,jk) = sa_tl(ji,jj,jk)
217                  END DO
218               END DO
219            END DO
220         ENDIF
221         !
222      ENDIF
223
224#if defined key_agrif
225      ! Update tracer at AGRIF zoom boundaries
226      error " Agrif not in tangent yet"
227#endif     
228
229      !
230   END SUBROUTINE tra_nxt_tan
231   SUBROUTINE tra_nxt_adj( kt )
232      !!----------------------------------------------------------------------
233      !!                   ***  ROUTINE tranxt_adj  ***
234      !!
235      !! ** Purpose of the direct routine:
236      !!             Apply the boundary condition on the after temperature 
237      !!             and salinity fields, achieved the time stepping by adding
238      !!             the Asselin filter on now fields and swapping the fields.
239      !!
240      !! ** Method  :   At this stage of the computation, ta and sa are the
241      !!             after temperature and salinity as the time stepping has
242      !!             been performed in trazdf_imp or trazdf_exp module.
243      !!
244      !!              - Apply lateral boundary conditions on (ta,sa)
245      !!             at the local domain   boundaries through lbc_lnk call,
246      !!             at the radiative open boundaries (lk_obc=T),
247      !!             at the relaxed   open boundaries (lk_bdy=T), and
248      !!             at the AGRIF zoom     boundaries (lk_agrif=T)
249      !!
250      !!              - Apply the Asselin time filter on now fields,
251      !!             save in (ta,sa) an average over the three time levels
252      !!             which will be used to compute rdn and thus the semi-implicit
253      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T), and
254      !!             swap tracer fields to prepare the next time_step.
255      !!                This can be summurized for tempearture as:
256      !!                    zt = (ta+2tn+tb)/4       ln_dynhpg_imp = T
257      !!                    zt = 0                   otherwise
258      !!                    tb = tn + atfp*[ tb - 2 tn + ta ]
259      !!                    tn = ta
260      !!                    ta = zt        (NB: reset to 0 after eos_bn2 call)
261      !!
262      !! ** Action  : - update (tb,sb) and (tn,sn)
263      !!              - (ta,sa) time averaged (t,s)      (ln_dynhpg_imp = T)
264      !!----------------------------------------------------------------------
265      !!
266      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
267      !!
268      !!
269      INTEGER  ::   ji, jj, jk    ! dummy loop indices
270      REAL(wp) ::   ztad, zsad    ! temporary scalars
271      REAL(wp) ::   zfact         ! temporary scalar
272
273      !!----------------------------------------------------------------------
274      IF( kt == nitend ) THEN
275         IF(lwp) WRITE(numout,*)
276         IF(lwp) WRITE(numout,*) 'tra_nxt_adj : achieve the time stepping by Asselin filter and array swap'
277         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
278      ENDIF
279
280#if defined key_agrif
281      ! Update tracer at AGRIF zoom boundaries
282      error " Agrif not in adjoint yet"
283#endif     
284      ! Asselin time filter and swap of arrays
285      !                                             
286      IF( neuler == 0 .AND. kt == nit000 ) THEN      ! Euler 1st time step : swap only
287         DO jk = 1, jpkm1
288            ta_ad(:,:,jk) = ta_ad(:,:,jk) + tn_ad(:,:,jk)
289            sa_ad(:,:,jk) = sa_ad(:,:,jk) + sn_ad(:,:,jk)
290
291            tn_ad(:,:,jk) = tb_ad(:,:,jk)
292            sn_ad(:,:,jk) = sb_ad(:,:,jk)
293            tb_ad(:,:,jk) = 0.0_wp
294            sb_ad(:,:,jk) = 0.0_wp
295         END DO
296         !                                           
297      ELSE                                           ! Leap-Frog : filter + swap
298         !                                           
299         IF( ln_dynhpg_imp ) THEN                         ! semi-implicite hpg case
300            DO jk = 1, jpkm1                              ! (save the averaged of the 3 time steps
301               DO jj = 1, jpj                             !                   in the after fields)
302                  DO ji = 1, jpi
303                     ztad = ta_ad(ji,jj,jk)
304                     zsad = sa_ad(ji,jj,jk)
305
306                     ta_ad(ji,jj,jk) =  tn_ad(ji,jj,jk) + tb_ad(ji,jj,jk) * atfp
307                     tn_ad(ji,jj,jk) =                    tb_ad(ji,jj,jk) * atfp1
308                     tb_ad(ji,jj,jk) =                    tb_ad(ji,jj,jk) * atfp
309
310                     sa_ad(ji,jj,jk) =  sn_ad(ji,jj,jk) + sb_ad(ji,jj,jk) * atfp 
311                     sn_ad(ji,jj,jk) =                    sb_ad(ji,jj,jk) * atfp1 
312                     sb_ad(ji,jj,jk) =                    sb_ad(ji,jj,jk) * atfp
313
314                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + ztad * 0.25_wp
315                     tn_ad(ji,jj,jk) = tn_ad(ji,jj,jk) + ztad * 0.5_wp
316                     tb_ad(ji,jj,jk) = tb_ad(ji,jj,jk) + ztad * 0.25_wp
317
318                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + zsad * 0.25_wp
319                     sn_ad(ji,jj,jk) = sn_ad(ji,jj,jk) + zsad * 0.5_wp
320                     sb_ad(ji,jj,jk) = sb_ad(ji,jj,jk) + zsad * 0.25_wp
321
322                  END DO
323               END DO
324            END DO
325         ELSE                                            ! explicit hpg case
326            DO jk = 1, jpkm1
327               DO jj = 1, jpj
328                  DO ji = 1, jpi
329                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + tn_ad(ji,jj,jk)
330                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + sn_ad(ji,jj,jk) 
331
332                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + atfp  * tb_ad(ji,jj,jk)
333                     tn_ad(ji,jj,jk) =                   atfp1 * tb_ad(ji,jj,jk)
334                     tb_ad(ji,jj,jk) =                   atfp  * tb_ad(ji,jj,jk)
335
336                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + atfp  * sb_ad(ji,jj,jk)
337                     sn_ad(ji,jj,jk) =                   atfp1 * sb_ad(ji,jj,jk)
338                     sb_ad(ji,jj,jk) =                   atfp  * sb_ad(ji,jj,jk)
339                  END DO
340               END DO
341            END DO
342         ENDIF
343         !
344      ENDIF
345#if defined key_agrif
346      error "agrif not available in tangent yet"
347#endif
348#if defined key_bdy
349      error "bdy not available in tangent yet"
350#endif
351#if defined key_obc
352      CALL obc_tra_adj( kt )               ! OBC open boundaries
353#endif
354      ! Update after tracer on domain lateral boundaries
355      !
356      CALL lbc_lnk_adj( ta_ad, 'T', 1. )      ! local domain boundaries  (T-point, unchanged sign)
357      CALL lbc_lnk_adj( sa_ad, 'T', 1. )
358      !
359   END SUBROUTINE tra_nxt_adj
360
361   SUBROUTINE tra_nxt_adj_tst( kumadt ) 
362      !!-----------------------------------------------------------------------
363      !!
364      !!          ***  ROUTINE tra_nxt_adj_tst : TEST OF tra_nxt_adj  ***
365      !!
366      !! ** Purpose : Test the adjoint routine.
367      !!
368      !! ** Method  : Verify the scalar product
369      !!           
370      !!                 ( L dx )^T W dy  =  dx^T L^T W dy
371      !!
372      !!              where  L   = tangent routine
373      !!                     L^T = adjoint routine
374      !!                     W   = diagonal matrix of scale factors
375      !!                     dx  = input perturbation (random field)
376      !!                     dy  = L dx
377      !!
378       !! History :
379      !!        ! 08-08 (A. Vidard)
380      !!-----------------------------------------------------------------------
381      !! * Modules used
382
383      !! * Arguments
384      INTEGER, INTENT(IN) :: &
385         & kumadt             ! Output unit
386 
387      INTEGER :: &
388         & ji,    &        ! dummy loop indices
389         & jj,    &       
390         & jk     
391      INTEGER, DIMENSION(jpi,jpj) :: &
392         & iseed_2d        ! 2D seed for the random number generator
393
394      !! * Local declarations
395      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: &
396         & zsb_tlin,     &! Tangent input : before salinity
397         & ztb_tlin,     &! Tangent input : before temperature
398         & zsa_tlin,     &! Tangent input : after salinity
399         & zta_tlin,     &! Tangent input : after temperature
400         & zsn_tlin,     &! Tangent input : now salinity
401         & ztn_tlin,     &! Tangent input : now temperature
402         & zsb_tlout,    &! Tangent output: before salinity
403         & ztb_tlout,    &! Tangent output: before temperature
404         & zsa_tlout,    &! Tangent output: after salinity
405         & zta_tlout,    &! Tangent output: after temperature
406         & zsn_tlout,    &! Tangent output: now salinity
407         & ztn_tlout,    &! Tangent output: now temperature
408         & zsb_adin,     &! Adjoint input : before salinity
409         & ztb_adin,     &! Adjoint input : before temperature
410         & zsa_adin,     &! Adjoint input : after salinity
411         & zta_adin,     &! Adjoint input : after temperature
412         & zsn_adin,     &! Adjoint input : now salinity
413         & ztn_adin,     &! Adjoint input : now temperature
414         & zsb_adout,    &! Adjoint output: before salinity
415         & ztb_adout,    &! Adjoint output: before temperature
416         & zsa_adout,    &! Adjoint output: after salinity
417         & zta_adout,    &! Adjoint output: after temperature
418         & zsn_adout,    &! Adjoint output: now salinity
419         & ztn_adout,    &! Adjoint output: now temperature
420         & zr             ! 3D field
421         
422      REAL(KIND=wp) ::       &
423         & zsp1,             & ! scalar product involving the tangent routine
424         & zsp1_1,           & ! scalar product involving the tangent routine
425         & zsp1_2,           & ! scalar product involving the tangent routine
426         & zsp1_3,           & ! scalar product involving the tangent routine
427         & zsp1_4,           & ! scalar product involving the tangent routine
428         & zsp1_5,           & ! scalar product involving the tangent routine
429         & zsp1_6,           & ! scalar product involving the tangent routine
430         & zsp2,             & ! scalar product involving the adjoint routine
431         & zsp2_1,           & ! scalar product involving the adjoint routine
432         & zsp2_2,           & ! scalar product involving the adjoint routine
433         & zsp2_3,           & ! scalar product involving the adjoint routine
434         & zsp2_4,           & ! scalar product involving the adjoint routine
435         & zsp2_5,           & ! scalar product involving the adjoint routine
436         & zsp2_6              ! scalar product involving the adjoint routine
437      CHARACTER(LEN=14) ::   &
438         & cl_name
439
440      ALLOCATE( & 
441         & zsb_tlin(jpi,jpj,jpk),    &! Tangent input : before salinity
442         & ztb_tlin(jpi,jpj,jpk),    &! Tangent input : before temperature
443         & zsa_tlin(jpi,jpj,jpk),    &! Tangent input : after salinity
444         & zta_tlin(jpi,jpj,jpk),    &! Tangent input : after temperature
445         & zsn_tlin(jpi,jpj,jpk),    &! Tangent input : now salinity
446         & ztn_tlin(jpi,jpj,jpk),    &! Tangent input : now temperature
447         & zsb_tlout(jpi,jpj,jpk),   &! Tangent output: before salinity
448         & ztb_tlout(jpi,jpj,jpk),   &! Tangent output: before temperature
449         & zsa_tlout(jpi,jpj,jpk),   &! Tangent output: after salinity
450         & zta_tlout(jpi,jpj,jpk),   &! Tangent output: after temperature
451         & zsn_tlout(jpi,jpj,jpk),   &! Tangent output: now salinity
452         & ztn_tlout(jpi,jpj,jpk),   &! Tangent output: now temperature
453         & zsb_adin(jpi,jpj,jpk),    &! Adjoint input : before salinity
454         & ztb_adin(jpi,jpj,jpk),    &! Adjoint input : before temperature
455         & zsa_adin(jpi,jpj,jpk),    &! Adjoint input : after salinity
456         & zta_adin(jpi,jpj,jpk),    &! Adjoint input : after temperature
457         & zsn_adin(jpi,jpj,jpk),    &! Adjoint input : now salinity
458         & ztn_adin(jpi,jpj,jpk),    &! Adjoint input : now temperature
459         & zsb_adout(jpi,jpj,jpk),   &! Adjoint output: before salinity
460         & ztb_adout(jpi,jpj,jpk),   &! Adjoint output: before temperature
461         & zsa_adout(jpi,jpj,jpk),   &! Adjoint output: after salinity
462         & zta_adout(jpi,jpj,jpk),   &! Adjoint output: after temperature
463         & zsn_adout(jpi,jpj,jpk),   &! Adjoint output: now salinity
464         & ztn_adout(jpi,jpj,jpk),   &! Adjoint output: now temperature
465         & zr       (jpi,jpj,jpk)    &! 3D field
466         & )
467
468
469      !==================================================================
470      ! 1) dx = ( tb_tl, tn_tl, ta_tl,       dy = ( tb_tl, tn_tl, ta_tl,
471      !           sb_tl, sn_tl, sa_tl  ) and        sb_tl, sn_tl, sa_tl )
472      !==================================================================
473
474      !--------------------------------------------------------------------
475      ! Reset the tangent and adjoint variables
476      !--------------------------------------------------------------------
477      sb_tl(:,:,:) = 0.0_wp
478      tb_tl(:,:,:) = 0.0_wp
479      sa_tl(:,:,:) = 0.0_wp
480      ta_tl(:,:,:) = 0.0_wp
481      sn_tl(:,:,:) = 0.0_wp
482      tn_tl(:,:,:) = 0.0_wp
483      sb_ad(:,:,:) = 0.0_wp
484      tb_ad(:,:,:) = 0.0_wp
485      sa_ad(:,:,:) = 0.0_wp
486      ta_ad(:,:,:) = 0.0_wp
487      sn_ad(:,:,:) = 0.0_wp
488      tn_ad(:,:,:) = 0.0_wp
489      zsb_tlin(:,:,:) = 0.0_wp
490      ztb_tlin(:,:,:) = 0.0_wp
491      zsa_tlin(:,:,:) = 0.0_wp
492      zta_tlin(:,:,:) = 0.0_wp
493      zsn_tlin(:,:,:) = 0.0_wp
494      ztn_tlin(:,:,:) = 0.0_wp
495
496      DO jj = 1, jpj
497         DO ji = 1, jpi
498            iseed_2d(ji,jj) = - ( 785483 + &
499               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
500         END DO
501      END DO
502      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
503      DO jk = 1, jpk
504         DO jj = nldj, nlej
505            DO ji = nldi, nlei
506               zsb_tlin(ji,jj,jk) = zr(ji,jj,jk) 
507            END DO
508         END DO
509      END DO
510
511      DO jj = 1, jpj
512         DO ji = 1, jpi
513            iseed_2d(ji,jj) = - ( 358606 + &
514               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
515         END DO
516      END DO
517      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
518      DO jk = 1, jpk
519         DO jj = nldj, nlej
520            DO ji = nldi, nlei
521               ztb_tlin(ji,jj,jk) = zr(ji,jj,jk) 
522            END DO
523         END DO
524      END DO
525
526      DO jj = 1, jpj
527         DO ji = 1, jpi
528            iseed_2d(ji,jj) = - ( 596035 + &
529               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
530         END DO
531      END DO
532      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
533      DO jk = 1, jpk
534         DO jj = nldj, nlej
535            DO ji = nldi, nlei
536               zsa_tlin(ji,jj,jk) = zr(ji,jj,jk) 
537            END DO
538         END DO
539      END DO
540
541      DO jj = 1, jpj
542         DO ji = 1, jpi
543            iseed_2d(ji,jj) = - ( 523432 + &
544               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
545         END DO
546      END DO
547      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
548      DO jk = 1, jpk
549         DO jj = nldj, nlej
550            DO ji = nldi, nlei
551               zta_tlin(ji,jj,jk) = zr(ji,jj,jk) 
552            END DO
553         END DO
554      END DO
555
556      DO jj = 1, jpj
557         DO ji = 1, jpi
558            iseed_2d(ji,jj) = - ( 263957 + &
559               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
560         END DO
561      END DO
562      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
563      DO jk = 1, jpk
564         DO jj = nldj, nlej
565            DO ji = nldi, nlei
566               zsn_tlin(ji,jj,jk) = zr(ji,jj,jk) 
567            END DO
568         END DO
569      END DO
570
571      DO jj = 1, jpj
572         DO ji = 1, jpi
573            iseed_2d(ji,jj) = - ( 459031 + &
574               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
575         END DO
576      END DO
577      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
578      DO jk = 1, jpk
579         DO jj = nldj, nlej
580            DO ji = nldi, nlei
581               ztn_tlin(ji,jj,jk) = zr(ji,jj,jk) 
582            END DO
583         END DO
584      END DO
585
586      sb_tl(:,:,:) = zsb_tlin(:,:,:) 
587      tb_tl(:,:,:) = ztb_tlin(:,:,:) 
588      sa_tl(:,:,:) = zsa_tlin(:,:,:) 
589      ta_tl(:,:,:) = zta_tlin(:,:,:)
590      sn_tl(:,:,:) = zsn_tlin(:,:,:) 
591      tn_tl(:,:,:) = ztn_tlin(:,:,:)
592
593      CALL tra_nxt_tan( nit000 + 1 )
594
595      zsa_tlout(:,:,:) = sa_tl(:,:,:) 
596      zta_tlout(:,:,:) = ta_tl(:,:,:)
597      zsb_tlout(:,:,:) = sb_tl(:,:,:) 
598      ztb_tlout(:,:,:) = tb_tl(:,:,:)
599      zsn_tlout(:,:,:) = sn_tl(:,:,:) 
600      ztn_tlout(:,:,:) = tn_tl(:,:,:)
601
602      !--------------------------------------------------------------------
603      ! Initialize the adjoint variables: dy^* = W dy
604      !--------------------------------------------------------------------
605
606      DO jk = 1, jpk
607        DO jj = nldj, nlej
608           DO ji = nldi, nlei
609              zsa_adin(ji,jj,jk)     = zsa_tlout(ji,jj,jk) &
610                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
611                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
612              zta_adin(ji,jj,jk)     = zta_tlout(ji,jj,jk) &
613                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
614                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
615              zsb_adin(ji,jj,jk)     = zsb_tlout(ji,jj,jk) &
616                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
617                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
618              ztb_adin(ji,jj,jk)     = ztb_tlout(ji,jj,jk) &
619                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
620                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
621              zsn_adin(ji,jj,jk)     = zsn_tlout(ji,jj,jk) &
622                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
623                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
624              ztn_adin(ji,jj,jk)     = ztn_tlout(ji,jj,jk) &
625                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
626                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
627            END DO
628         END DO
629      END DO
630
631      !--------------------------------------------------------------------
632      ! Compute the scalar product: ( L dx )^T W dy
633      !--------------------------------------------------------------------
634
635      zsp1_1 = DOT_PRODUCT( zsa_tlout    , zsa_adin     )
636      zsp1_2 = DOT_PRODUCT( zta_tlout    , zta_adin     )
637      zsp1_3 = DOT_PRODUCT( zsb_tlout    , zsb_adin     )
638      zsp1_4 = DOT_PRODUCT( ztb_tlout    , ztb_adin     )
639      zsp1_5 = DOT_PRODUCT( zsn_tlout    , zsn_adin     )
640      zsp1_6 = DOT_PRODUCT( ztn_tlout    , ztn_adin     )
641      zsp1   = zsp1_1 + zsp1_2 + zsp1_3 + zsp1_4 + zsp1_5 + zsp1_6
642
643      !--------------------------------------------------------------------
644      ! Call the adjoint routine: dx^* = L^T dy^*
645      !--------------------------------------------------------------------
646
647      sa_ad(:,:,:)     = zsa_adin(:,:,:)
648      ta_ad(:,:,:)     = zta_adin(:,:,:)
649      sb_ad(:,:,:)     = zsb_adin(:,:,:)
650      tb_ad(:,:,:)     = ztb_adin(:,:,:)
651      sn_ad(:,:,:)     = zsn_adin(:,:,:)
652      tn_ad(:,:,:)     = ztn_adin(:,:,:)
653
654      CALL tra_nxt_adj ( nit000 + 1 )
655
656      zsb_adout(:,:,:) = sb_ad(:,:,:)
657      ztb_adout(:,:,:) = tb_ad(:,:,:)
658      zsa_adout(:,:,:) = sa_ad(:,:,:)
659      zta_adout(:,:,:) = ta_ad(:,:,:)
660      zsn_adout(:,:,:) = sn_ad(:,:,:)
661      ztn_adout(:,:,:) = tn_ad(:,:,:)
662
663      !--------------------------------------------------------------------
664      ! Compute the scalar product: dx^T L^T W dy
665      !--------------------------------------------------------------------
666
667      zsp2_1 = DOT_PRODUCT( zsb_tlin  , zsb_adout   )
668      zsp2_2 = DOT_PRODUCT( ztb_tlin  , ztb_adout   )
669      zsp2_3 = DOT_PRODUCT( zsa_tlin  , zsa_adout   )
670      zsp2_4 = DOT_PRODUCT( zta_tlin  , zta_adout   )
671      zsp2_5 = DOT_PRODUCT( zsn_tlin  , zsn_adout   )
672      zsp2_6 = DOT_PRODUCT( ztn_tlin  , ztn_adout   )
673
674      zsp2 = zsp2_1 + zsp2_2 + zsp2_3 + zsp2_4 + zsp2_5 + zsp2_6
675
676      ! Compare the scalar products
677
678      ! 14 char:'12345678901234'
679      cl_name = 'tra_nxt_adj   '
680      CALL prntst_adj( cl_name, kumadt, zsp1, zsp2 )
681
682
683      DEALLOCATE( & 
684         & zsb_tlin,     &
685         & ztb_tlin,     &
686         & zsa_tlin,     &
687         & zta_tlin,     &
688         & zsn_tlin,     &
689         & ztn_tlin,     &
690         & zsb_tlout,    &
691         & ztb_tlout,    &
692         & zsa_tlout,    &
693         & zta_tlout,    &
694         & zsn_tlout,    &
695         & ztn_tlout,    &
696         & zsb_adin,     &
697         & ztb_adin,     &
698         & zsa_adin,     &
699         & zta_adin,     &
700         & zsn_adin,     &
701         & ztn_adin,     &
702         & zsb_adout,    &
703         & ztb_adout,    &
704         & zsa_adout,    &
705         & zta_adout,    &
706         & zsn_adout,    &
707         & ztn_adout,    &
708         & zr            &
709         & )
710
711     END SUBROUTINE tra_nxt_adj_tst
712   !!======================================================================
713#endif
714END MODULE tranxt_tam
Note: See TracBrowser for help on using the repository browser.