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 @ 1885

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

add TAM sources

  • Property svn:executable set to *
File size: 30.0 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
232   SUBROUTINE tra_nxt_adj( kt )
233      !!----------------------------------------------------------------------
234      !!                   ***  ROUTINE tranxt_adj  ***
235      !!
236      !! ** Purpose of the direct routine:
237      !!             Apply the boundary condition on the after temperature 
238      !!             and salinity fields, achieved the time stepping by adding
239      !!             the Asselin filter on now fields and swapping the fields.
240      !!
241      !! ** Method  :   At this stage of the computation, ta and sa are the
242      !!             after temperature and salinity as the time stepping has
243      !!             been performed in trazdf_imp or trazdf_exp module.
244      !!
245      !!              - Apply lateral boundary conditions on (ta,sa)
246      !!             at the local domain   boundaries through lbc_lnk call,
247      !!             at the radiative open boundaries (lk_obc=T),
248      !!             at the relaxed   open boundaries (lk_bdy=T), and
249      !!             at the AGRIF zoom     boundaries (lk_agrif=T)
250      !!
251      !!              - Apply the Asselin time filter on now fields,
252      !!             save in (ta,sa) an average over the three time levels
253      !!             which will be used to compute rdn and thus the semi-implicit
254      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T), and
255      !!             swap tracer fields to prepare the next time_step.
256      !!                This can be summurized for tempearture as:
257      !!                    zt = (ta+2tn+tb)/4       ln_dynhpg_imp = T
258      !!                    zt = 0                   otherwise
259      !!                    tb = tn + atfp*[ tb - 2 tn + ta ]
260      !!                    tn = ta
261      !!                    ta = zt        (NB: reset to 0 after eos_bn2 call)
262      !!
263      !! ** Action  : - update (tb,sb) and (tn,sn)
264      !!              - (ta,sa) time averaged (t,s)      (ln_dynhpg_imp = T)
265      !!----------------------------------------------------------------------
266      !!
267      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
268      !!
269      !!
270      INTEGER  ::   ji, jj, jk    ! dummy loop indices
271      REAL(wp) ::   ztad, zsad    ! temporary scalars
272      REAL(wp) ::   zfact         ! temporary scalar
273
274      !!----------------------------------------------------------------------
275      IF( kt == nitend ) THEN
276         IF(lwp) WRITE(numout,*)
277         IF(lwp) WRITE(numout,*) 'tra_nxt_adj : achieve the time stepping by Asselin filter and array swap'
278         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
279      ENDIF
280
281#if defined key_agrif
282      ! Update tracer at AGRIF zoom boundaries
283      error " Agrif not in adjoint yet"
284#endif     
285      ! Asselin time filter and swap of arrays
286      !                                             
287      IF( neuler == 0 .AND. kt == nit000 ) THEN      ! Euler 1st time step : swap only
288         DO jk = 1, jpkm1
289            ta_ad(:,:,jk) = ta_ad(:,:,jk) + tn_ad(:,:,jk)
290            sa_ad(:,:,jk) = sa_ad(:,:,jk) + sn_ad(:,:,jk)
291
292            tn_ad(:,:,jk) = tb_ad(:,:,jk)
293            sn_ad(:,:,jk) = sb_ad(:,:,jk)
294            tb_ad(:,:,jk) = 0.0_wp
295            sb_ad(:,:,jk) = 0.0_wp
296         END DO
297         !                                           
298      ELSE                                           ! Leap-Frog : filter + swap
299         !                                           
300         IF( ln_dynhpg_imp ) THEN                         ! semi-implicite hpg case
301            DO jk = 1, jpkm1                              ! (save the averaged of the 3 time steps
302               DO jj = 1, jpj                             !                   in the after fields)
303                  DO ji = 1, jpi
304                     ztad = ta_ad(ji,jj,jk)
305                     zsad = sa_ad(ji,jj,jk)
306
307                     ta_ad(ji,jj,jk) =  tn_ad(ji,jj,jk) + tb_ad(ji,jj,jk) * atfp
308                     tn_ad(ji,jj,jk) =                    tb_ad(ji,jj,jk) * atfp1
309                     tb_ad(ji,jj,jk) =                    tb_ad(ji,jj,jk) * atfp
310
311                     sa_ad(ji,jj,jk) =  sn_ad(ji,jj,jk) + sb_ad(ji,jj,jk) * atfp 
312                     sn_ad(ji,jj,jk) =                    sb_ad(ji,jj,jk) * atfp1 
313                     sb_ad(ji,jj,jk) =                    sb_ad(ji,jj,jk) * atfp
314
315                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + ztad * 0.25_wp
316                     tn_ad(ji,jj,jk) = tn_ad(ji,jj,jk) + ztad * 0.5_wp
317                     tb_ad(ji,jj,jk) = tb_ad(ji,jj,jk) + ztad * 0.25_wp
318
319                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + zsad * 0.25_wp
320                     sn_ad(ji,jj,jk) = sn_ad(ji,jj,jk) + zsad * 0.5_wp
321                     sb_ad(ji,jj,jk) = sb_ad(ji,jj,jk) + zsad * 0.25_wp
322
323                  END DO
324               END DO
325            END DO
326         ELSE                                            ! explicit hpg case
327            DO jk = 1, jpkm1
328               DO jj = 1, jpj
329                  DO ji = 1, jpi
330                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + tn_ad(ji,jj,jk)
331                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + sn_ad(ji,jj,jk) 
332
333                     ta_ad(ji,jj,jk) = ta_ad(ji,jj,jk) + atfp  * tb_ad(ji,jj,jk)
334                     tn_ad(ji,jj,jk) =                   atfp1 * tb_ad(ji,jj,jk)
335                     tb_ad(ji,jj,jk) =                   atfp  * tb_ad(ji,jj,jk)
336
337                     sa_ad(ji,jj,jk) = sa_ad(ji,jj,jk) + atfp  * sb_ad(ji,jj,jk)
338                     sn_ad(ji,jj,jk) =                   atfp1 * sb_ad(ji,jj,jk)
339                     sb_ad(ji,jj,jk) =                   atfp  * sb_ad(ji,jj,jk)
340                  END DO
341               END DO
342            END DO
343         ENDIF
344         !
345      ENDIF
346#if defined key_agrif
347      error "agrif not available in tangent yet"
348#endif
349#if defined key_bdy
350      error "bdy not available in tangent yet"
351#endif
352#if defined key_obc
353      CALL obc_tra_adj( kt )               ! OBC open boundaries
354#endif
355      ! Update after tracer on domain lateral boundaries
356      !
357      CALL lbc_lnk_adj( ta_ad, 'T', 1. )      ! local domain boundaries  (T-point, unchanged sign)
358      CALL lbc_lnk_adj( sa_ad, 'T', 1. )
359      !
360   END SUBROUTINE tra_nxt_adj
361
362   SUBROUTINE tra_nxt_adj_tst( kumadt ) 
363      !!-----------------------------------------------------------------------
364      !!
365      !!          ***  ROUTINE tra_nxt_adj_tst : TEST OF tra_nxt_adj  ***
366      !!
367      !! ** Purpose : Test the adjoint routine.
368      !!
369      !! ** Method  : Verify the scalar product
370      !!           
371      !!                 ( L dx )^T W dy  =  dx^T L^T W dy
372      !!
373      !!              where  L   = tangent routine
374      !!                     L^T = adjoint routine
375      !!                     W   = diagonal matrix of scale factors
376      !!                     dx  = input perturbation (random field)
377      !!                     dy  = L dx
378      !!
379       !! History :
380      !!        ! 08-08 (A. Vidard)
381      !!-----------------------------------------------------------------------
382      !! * Modules used
383
384      !! * Arguments
385      INTEGER, INTENT(IN) :: &
386         & kumadt             ! Output unit
387 
388      INTEGER :: &
389         & ji,    &        ! dummy loop indices
390         & jj,    &       
391         & jk     
392      INTEGER, DIMENSION(jpi,jpj) :: &
393         & iseed_2d        ! 2D seed for the random number generator
394
395      !! * Local declarations
396      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: &
397         & zsb_tlin,     &! Tangent input : before salinity
398         & ztb_tlin,     &! Tangent input : before temperature
399         & zsa_tlin,     &! Tangent input : after salinity
400         & zta_tlin,     &! Tangent input : after temperature
401         & zsn_tlin,     &! Tangent input : now salinity
402         & ztn_tlin,     &! Tangent input : now temperature
403         & zsb_tlout,    &! Tangent output: before salinity
404         & ztb_tlout,    &! Tangent output: before temperature
405         & zsa_tlout,    &! Tangent output: after salinity
406         & zta_tlout,    &! Tangent output: after temperature
407         & zsn_tlout,    &! Tangent output: now salinity
408         & ztn_tlout,    &! Tangent output: now temperature
409         & zsb_adin,     &! Adjoint input : before salinity
410         & ztb_adin,     &! Adjoint input : before temperature
411         & zsa_adin,     &! Adjoint input : after salinity
412         & zta_adin,     &! Adjoint input : after temperature
413         & zsn_adin,     &! Adjoint input : now salinity
414         & ztn_adin,     &! Adjoint input : now temperature
415         & zsb_adout,    &! Adjoint output: before salinity
416         & ztb_adout,    &! Adjoint output: before temperature
417         & zsa_adout,    &! Adjoint output: after salinity
418         & zta_adout,    &! Adjoint output: after temperature
419         & zsn_adout,    &! Adjoint output: now salinity
420         & ztn_adout,    &! Adjoint output: now temperature
421         & zr             ! 3D field
422         
423      REAL(KIND=wp) ::       &
424         & zsp1,             & ! scalar product involving the tangent routine
425         & zsp1_1,           & ! scalar product involving the tangent routine
426         & zsp1_2,           & ! scalar product involving the tangent routine
427         & zsp1_3,           & ! scalar product involving the tangent routine
428         & zsp1_4,           & ! scalar product involving the tangent routine
429         & zsp1_5,           & ! scalar product involving the tangent routine
430         & zsp1_6,           & ! scalar product involving the tangent routine
431         & zsp2,             & ! scalar product involving the adjoint routine
432         & zsp2_1,           & ! scalar product involving the adjoint routine
433         & zsp2_2,           & ! scalar product involving the adjoint routine
434         & zsp2_3,           & ! scalar product involving the adjoint routine
435         & zsp2_4,           & ! scalar product involving the adjoint routine
436         & zsp2_5,           & ! scalar product involving the adjoint routine
437         & zsp2_6              ! scalar product involving the adjoint routine
438      CHARACTER(LEN=14) ::   &
439         & cl_name
440
441      ALLOCATE( & 
442         & zsb_tlin(jpi,jpj,jpk),    &! Tangent input : before salinity
443         & ztb_tlin(jpi,jpj,jpk),    &! Tangent input : before temperature
444         & zsa_tlin(jpi,jpj,jpk),    &! Tangent input : after salinity
445         & zta_tlin(jpi,jpj,jpk),    &! Tangent input : after temperature
446         & zsn_tlin(jpi,jpj,jpk),    &! Tangent input : now salinity
447         & ztn_tlin(jpi,jpj,jpk),    &! Tangent input : now temperature
448         & zsb_tlout(jpi,jpj,jpk),   &! Tangent output: before salinity
449         & ztb_tlout(jpi,jpj,jpk),   &! Tangent output: before temperature
450         & zsa_tlout(jpi,jpj,jpk),   &! Tangent output: after salinity
451         & zta_tlout(jpi,jpj,jpk),   &! Tangent output: after temperature
452         & zsn_tlout(jpi,jpj,jpk),   &! Tangent output: now salinity
453         & ztn_tlout(jpi,jpj,jpk),   &! Tangent output: now temperature
454         & zsb_adin(jpi,jpj,jpk),    &! Adjoint input : before salinity
455         & ztb_adin(jpi,jpj,jpk),    &! Adjoint input : before temperature
456         & zsa_adin(jpi,jpj,jpk),    &! Adjoint input : after salinity
457         & zta_adin(jpi,jpj,jpk),    &! Adjoint input : after temperature
458         & zsn_adin(jpi,jpj,jpk),    &! Adjoint input : now salinity
459         & ztn_adin(jpi,jpj,jpk),    &! Adjoint input : now temperature
460         & zsb_adout(jpi,jpj,jpk),   &! Adjoint output: before salinity
461         & ztb_adout(jpi,jpj,jpk),   &! Adjoint output: before temperature
462         & zsa_adout(jpi,jpj,jpk),   &! Adjoint output: after salinity
463         & zta_adout(jpi,jpj,jpk),   &! Adjoint output: after temperature
464         & zsn_adout(jpi,jpj,jpk),   &! Adjoint output: now salinity
465         & ztn_adout(jpi,jpj,jpk),   &! Adjoint output: now temperature
466         & zr       (jpi,jpj,jpk)    &! 3D field
467         & )
468
469
470      !==================================================================
471      ! 1) dx = ( tb_tl, tn_tl, ta_tl,       dy = ( tb_tl, tn_tl, ta_tl,
472      !           sb_tl, sn_tl, sa_tl  ) and        sb_tl, sn_tl, sa_tl )
473      !==================================================================
474
475      !--------------------------------------------------------------------
476      ! Reset the tangent and adjoint variables
477      !--------------------------------------------------------------------
478      zsa_tlin(:,:,:) = 0.0_wp 
479      zta_tlin(:,:,:) = 0.0_wp 
480      zsb_tlin(:,:,:) = 0.0_wp 
481      ztb_tlin(:,:,:) = 0.0_wp 
482      zsn_tlin(:,:,:) = 0.0_wp 
483      ztn_tlin(:,:,:) = 0.0_wp 
484      zsa_adin(:,:,:) = 0.0_wp 
485      zta_adin(:,:,:) = 0.0_wp 
486      zsb_adin(:,:,:) = 0.0_wp 
487      ztb_adin(:,:,:) = 0.0_wp 
488      zsn_adin(:,:,:) = 0.0_wp 
489      ztn_adin(:,:,:) = 0.0_wp 
490      sb_tl(:,:,:) = 0.0_wp
491      tb_tl(:,:,:) = 0.0_wp
492      sa_tl(:,:,:) = 0.0_wp
493      ta_tl(:,:,:) = 0.0_wp
494      sn_tl(:,:,:) = 0.0_wp
495      tn_tl(:,:,:) = 0.0_wp
496      sb_ad(:,:,:) = 0.0_wp
497      tb_ad(:,:,:) = 0.0_wp
498      sa_ad(:,:,:) = 0.0_wp
499      ta_ad(:,:,:) = 0.0_wp
500      sn_ad(:,:,:) = 0.0_wp
501      tn_ad(:,:,:) = 0.0_wp
502
503      DO jj = 1, jpj
504         DO ji = 1, jpi
505            iseed_2d(ji,jj) = - ( 785483 + &
506               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
507         END DO
508      END DO
509      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
510      DO jk = 1, jpk
511         DO jj = nldj, nlej
512            DO ji = nldi, nlei
513               zsb_tlin(ji,jj,jk) = zr(ji,jj,jk) 
514            END DO
515         END DO
516      END DO
517
518      DO jj = 1, jpj
519         DO ji = 1, jpi
520            iseed_2d(ji,jj) = - ( 358606 + &
521               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
522         END DO
523      END DO
524      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
525      DO jk = 1, jpk
526         DO jj = nldj, nlej
527            DO ji = nldi, nlei
528               ztb_tlin(ji,jj,jk) = zr(ji,jj,jk) 
529            END DO
530         END DO
531      END DO
532
533      DO jj = 1, jpj
534         DO ji = 1, jpi
535            iseed_2d(ji,jj) = - ( 596035 + &
536               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
537         END DO
538      END DO
539      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
540      DO jk = 1, jpk
541         DO jj = nldj, nlej
542            DO ji = nldi, nlei
543               zsa_tlin(ji,jj,jk) = zr(ji,jj,jk) 
544            END DO
545         END DO
546      END DO
547
548      DO jj = 1, jpj
549         DO ji = 1, jpi
550            iseed_2d(ji,jj) = - ( 523432 + &
551               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
552         END DO
553      END DO
554      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
555      DO jk = 1, jpk
556         DO jj = nldj, nlej
557            DO ji = nldi, nlei
558               zta_tlin(ji,jj,jk) = zr(ji,jj,jk) 
559            END DO
560         END DO
561      END DO
562
563      DO jj = 1, jpj
564         DO ji = 1, jpi
565            iseed_2d(ji,jj) = - ( 263957 + &
566               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
567         END DO
568      END DO
569      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stds )
570      DO jk = 1, jpk
571         DO jj = nldj, nlej
572            DO ji = nldi, nlei
573               zsn_tlin(ji,jj,jk) = zr(ji,jj,jk) 
574            END DO
575         END DO
576      END DO
577
578      DO jj = 1, jpj
579         DO ji = 1, jpi
580            iseed_2d(ji,jj) = - ( 459031 + &
581               &                  mig(ji) + ( mjg(jj) - 1 ) * jpiglo )
582         END DO
583      END DO
584      CALL grid_random( iseed_2d, zr, 'T', 0.0_wp, stdt )
585      DO jk = 1, jpk
586         DO jj = nldj, nlej
587            DO ji = nldi, nlei
588               ztn_tlin(ji,jj,jk) = zr(ji,jj,jk) 
589            END DO
590         END DO
591      END DO
592
593      sb_tl(:,:,:) = zsb_tlin(:,:,:) 
594      tb_tl(:,:,:) = ztb_tlin(:,:,:) 
595      sa_tl(:,:,:) = zsa_tlin(:,:,:) 
596      ta_tl(:,:,:) = zta_tlin(:,:,:)
597      sn_tl(:,:,:) = zsn_tlin(:,:,:) 
598      tn_tl(:,:,:) = ztn_tlin(:,:,:)
599
600      CALL tra_nxt_tan( nit000 + 1 )
601
602      zsa_tlout(:,:,:) = sa_tl(:,:,:) 
603      zta_tlout(:,:,:) = ta_tl(:,:,:)
604      zsb_tlout(:,:,:) = sb_tl(:,:,:) 
605      ztb_tlout(:,:,:) = tb_tl(:,:,:)
606      zsn_tlout(:,:,:) = sn_tl(:,:,:) 
607      ztn_tlout(:,:,:) = tn_tl(:,:,:)
608
609      !--------------------------------------------------------------------
610      ! Initialize the adjoint variables: dy^* = W dy
611      !--------------------------------------------------------------------
612
613      DO jk = 1, jpk
614        DO jj = nldj, nlej
615           DO ji = nldi, nlei
616              zsa_adin(ji,jj,jk)     = zsa_tlout(ji,jj,jk) &
617                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
618                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
619              zta_adin(ji,jj,jk)     = zta_tlout(ji,jj,jk) &
620                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
621                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
622              zsb_adin(ji,jj,jk)     = zsb_tlout(ji,jj,jk) &
623                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
624                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
625              ztb_adin(ji,jj,jk)     = ztb_tlout(ji,jj,jk) &
626                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
627                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
628              zsn_adin(ji,jj,jk)     = zsn_tlout(ji,jj,jk) &
629                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
630                 &                   * tmask(ji,jj,jk) * wesp_s(jk)
631              ztn_adin(ji,jj,jk)     = ztn_tlout(ji,jj,jk) &
632                 &                   * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) &
633                 &                   * tmask(ji,jj,jk) * wesp_t(jk)
634            END DO
635         END DO
636      END DO
637
638      !--------------------------------------------------------------------
639      ! Compute the scalar product: ( L dx )^T W dy
640      !--------------------------------------------------------------------
641
642      zsp1_1 = DOT_PRODUCT( zsa_tlout    , zsa_adin     )
643      zsp1_2 = DOT_PRODUCT( zta_tlout    , zta_adin     )
644      zsp1_3 = DOT_PRODUCT( zsb_tlout    , zsb_adin     )
645      zsp1_4 = DOT_PRODUCT( ztb_tlout    , ztb_adin     )
646      zsp1_5 = DOT_PRODUCT( zsn_tlout    , zsn_adin     )
647      zsp1_6 = DOT_PRODUCT( ztn_tlout    , ztn_adin     )
648      zsp1   = zsp1_1 + zsp1_2 + zsp1_3 + zsp1_4 + zsp1_5 + zsp1_6
649
650      !--------------------------------------------------------------------
651      ! Call the adjoint routine: dx^* = L^T dy^*
652      !--------------------------------------------------------------------
653
654      sa_ad(:,:,:)     = zsa_adin(:,:,:)
655      ta_ad(:,:,:)     = zta_adin(:,:,:)
656      sb_ad(:,:,:)     = zsb_adin(:,:,:)
657      tb_ad(:,:,:)     = ztb_adin(:,:,:)
658      sn_ad(:,:,:)     = zsn_adin(:,:,:)
659      tn_ad(:,:,:)     = ztn_adin(:,:,:)
660
661      CALL tra_nxt_adj ( nit000 + 1 )
662
663      zsb_adout(:,:,:) = sb_ad(:,:,:)
664      ztb_adout(:,:,:) = tb_ad(:,:,:)
665      zsa_adout(:,:,:) = sa_ad(:,:,:)
666      zta_adout(:,:,:) = ta_ad(:,:,:)
667      zsn_adout(:,:,:) = sn_ad(:,:,:)
668      ztn_adout(:,:,:) = tn_ad(:,:,:)
669
670      !--------------------------------------------------------------------
671      ! Compute the scalar product: dx^T L^T W dy
672      !--------------------------------------------------------------------
673
674      zsp2_1 = DOT_PRODUCT( zsb_tlin  , zsb_adout   )
675      zsp2_2 = DOT_PRODUCT( ztb_tlin  , ztb_adout   )
676      zsp2_3 = DOT_PRODUCT( zsa_tlin  , zsa_adout   )
677      zsp2_4 = DOT_PRODUCT( zta_tlin  , zta_adout   )
678      zsp2_5 = DOT_PRODUCT( zsn_tlin  , zsn_adout   )
679      zsp2_6 = DOT_PRODUCT( ztn_tlin  , ztn_adout   )
680
681      zsp2 = zsp2_1 + zsp2_2 + zsp2_3 + zsp2_4 + zsp2_5 + zsp2_6
682
683      ! Compare the scalar products
684
685      ! 14 char:'12345678901234'
686      cl_name = 'tra_nxt_adj   '
687      CALL prntst_adj( cl_name, kumadt, zsp1, zsp2 )
688
689
690      DEALLOCATE( & 
691         & zsb_tlin,     &
692         & ztb_tlin,     &
693         & zsa_tlin,     &
694         & zta_tlin,     &
695         & zsn_tlin,     &
696         & ztn_tlin,     &
697         & zsb_tlout,    &
698         & ztb_tlout,    &
699         & zsa_tlout,    &
700         & zta_tlout,    &
701         & zsn_tlout,    &
702         & ztn_tlout,    &
703         & zsb_adin,     &
704         & ztb_adin,     &
705         & zsa_adin,     &
706         & zta_adin,     &
707         & zsn_adin,     &
708         & ztn_adin,     &
709         & zsb_adout,    &
710         & ztb_adout,    &
711         & zsa_adout,    &
712         & zta_adout,    &
713         & zsn_adout,    &
714         & ztn_adout,    &
715         & zr            &
716         & )
717
718     END SUBROUTINE tra_nxt_adj_tst
719   !!======================================================================
720#endif
721END MODULE tranxt_tam
Note: See TracBrowser for help on using the repository browser.