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.
tamtst.F90 in branches/2012/dev_v3_4_STABLE_2012/NEMOGCM/NEMO/OPATAM_SRC – NEMO

source: branches/2012/dev_v3_4_STABLE_2012/NEMOGCM/NEMO/OPATAM_SRC/tamtst.F90 @ 4600

Last change on this file since 4600 was 3658, checked in by pabouttier, 12 years ago

Missing allocation for sbcssr_tam variables - see Ticket #1022

  • Property svn:executable set to *
File size: 10.9 KB
Line 
1MODULE tamtst
2#if defined key_tam
3   !!======================================================================
4   !!                       ***  MODULE tst ***
5   !! NEMOVAR : Testing routine
6   !!======================================================================
7
8   !!----------------------------------------------------------------------
9   !!----------------------------------------------------------------------
10   !!----------------------------------------------------------------------
11   !! * Modules used
12   USE par_oce
13   USE dom_oce
14   USE trabbc
15   USE traqsr
16   USE trabbl
17   USE in_out_manager      ! Input/output
18   USE iom                 ! netCDF output/input
19   USE paresp              ! Weights for energy-type scalar product
20   USE tamctl              ! TAM control
21   USE tradmp
22   USE sbcmod_tam          ! Tangent/adjoint of surface BCs
23   USE eosbn2_tam          ! Tangent/adjoint of eq. of state, Brunt-Vaisala
24   USE trasbc_tam          ! Tangent/adjoint of surface BCs application
25   USE traqsr_tam          ! Tangent/adjoint of penetrative solar radiation
26   USE trabbc_tam          ! Tangent/adjoint of bottom heat flux
27   USE trabbl_tam          ! Tangent/adjoint of bottom boundary layer
28   USE tradmp_tam          ! Tangent/adjoint of internal damping trends
29   USE traadv_tam          ! Tangent/adjoint of horizontal/vertical advection
30   USE cla_tam             ! Tangent/adjoint of cross land advection
31   USE traldf_tam          ! Tangent/adjoint of lateral mixing
32   USE trazdf_tam          ! Tangent/adjoint of vertical diffusion
33   USE tranxt_tam          ! Tangent/adjoint of tracers at next time step
34   USE zpshde_tam          ! Tangent/adjoint of horiz. derivs. for partial steps
35   USE divcur_tam          ! Tangent/adjoint of horiz. div. and rel. vorticity
36   USE dynadv_tam          ! Tangent/adjoint of horizontal/vertical advection
37   USE dynhpg_tam          ! Tangent/adjoint of horiz. pressure gradient
38   USE dynkeg_tam          ! Tangent/adjoint of kinetic energy gradient
39   USE dynldf_tam          ! Tangent/adjoint of lateral mixing
40   USE dynnxt_tam          ! Tangent/adjoint of dynamics at next time step
41   USE dynspg_tam          ! Tangent/adjoint of surface pressure gradient
42   USE dynvor_tam          ! Tangent/adjoint of relative and planetary vorticity
43   USE dynzdf_tam          ! Tangent/adjoint of vertical diffusion
44   USE dynbfr_tam          ! Tangent/adjoint of bottom friction
45   USE sshwzv_tam          ! Tangent/adjoint of vertical velocity
46   USE step_tam            ! manager of the adjoint ocean time stepping
47   USE trj_tam             ! reference trajectory
48   USE istate_tam
49   USE solsor_tam
50   USE closea_tam
51   USE zdfbfr_tam
52#if defined key_mpp_mpi
53   USE lbcnfd_tam
54#endif
55
56   IMPLICIT NONE
57
58   !! * Routine accessibility
59   PRIVATE
60
61   PUBLIC &
62      & tam_tst,  &        !: Scalar product test of the adjoint routines
63      & tam_tst_init, &    !: Reading of the namelist
64      & numadt             !: File unit number for adjoint test output
65
66   !! * Module variables
67   INTEGER :: &
68      & numadt             !: File unit number for adjoint test output
69
70CONTAINS
71
72   SUBROUTINE tam_tst
73      !!-----------------------------------------------------------------------
74      !!
75      !!                  ***  ROUTINE tam_tst  ***
76      !!
77      !! ** Purpose : Apply various tests (linearization, adjoint)
78      !!              on the NEMOTAM code.
79      !!
80      !! ** Method  :
81      !!
82      !! ** Action  :
83      !!
84      !! History :
85      !!        ! 2007-11 (A. Weaver) original (adjoint tests)
86      !!        ! 2009-08 (F. Vigilant) Add tangent tests
87      !!-----------------------------------------------------------------------
88      !! * Modules used
89
90      !! * Arguments
91
92      !! * Local declarations
93      CHARACTER (LEN=128) :: file_out
94
95      ! Open adjoint test output unit
96
97      IF (lwp) THEN
98
99         WRITE(file_out,FMT="('adjoint_test.output_',I4.4)") &
100            &   narea-1
101         CALL ctl_opn( numadt, file_out, 'UNKNOWN', 'FORMATTED',   &
102            &         'SEQUENTIAL', 1, numadt, .FALSE., 1 )
103
104         WRITE(numout,*) ' tstopt: Start testing adjoint operators ...'
105         WRITE(numout,*) ' ------'
106
107         WRITE(numout,*)
108         WRITE(numout,990) file_out
109990      FORMAT('          Output in file = ',A20)
110         WRITE(numout,*)
111
112         WRITE(numadt,*)
113         WRITE(numadt,997)
114         WRITE(numadt,998)
115         WRITE(numadt,999)
116997      FORMAT('  Routine (L)',2X,' ( L * dx )^T W dy ',2X, &
117            &   '     dx^T L^T W dy    ',2X,'Rel.',2X,      &
118            &   'Mach.',2X,'Status')
119998      FORMAT('             ',2X,'                   ',2X, &
120            &   '                      ',2X,'Err.',2x,      &
121            &   'Eps. ',2X,'      ')
122999      FORMAT('  -----------',2X,'-------------------',2X, &
123            &   '----------------------',2X,'----',2X,      &
124            &   '-----',2X,'------')
125         CALL FLUSH(numout)
126         CALL FLUSH(numadt)
127
128      ENDIF
129
130      ! Initialize energy weights
131      CALL par_esp
132
133      ! -----------------------------------------------------
134      ! Test the adjoint of the components of M (NEMOTAM)
135      ! -----------------------------------------------------
136
137      IF ( ln_swi_opatam == 0 ) THEN
138         !
139         ! *** initialize the reference trajectory
140         ! ------------
141          CALL trj_rea( nit000 - 1, 1 )
142         ! *** Tracers
143         ! -----------
144         ! *** Surface boundary conditions
145         ! ------------
146         CALL sbc_adj_tst( numadt )        ! surface boundary conditions
147         CALL tra_adv_adj_tst( numadt )    ! Horizontal and vertical advection
148         IF (lwp) WRITE(numadt,*)
149         IF ( ln_traqsr      )  THEN
150            CALL tra_qsr_adj_tst( numadt )    ! Penetrative solar radiation
151            IF (lwp) WRITE(numadt,*)
152         ENDIF
153         CALL tra_ldf_adj_tst( numadt )    ! Lateral mixing
154         IF (lwp) WRITE(numadt,*)
155         CALL eos_adj_tst( numadt )        ! In situ density
156         IF (lwp) WRITE(numadt,*)
157         IF ( ln_tradmp      )  THEN
158            CALL tra_dmp_adj_tst( numadt )    ! Internal damping trends
159            IF (lwp) WRITE(numadt,*)
160         ENDIF
161# if defined key_trabbl   ||   defined key_esopa
162         CALL tra_bbl_adj_tst( numadt )! Bottom boundary layer
163         IF (lwp) WRITE(numadt,*)
164# endif
165         IF ( nn_cla == 1     )  THEN
166            CALL cla_traadv_adj_tst( numadt )    ! Cross land advection (update hor. advection)
167            IF (lwp) WRITE(numadt,*)
168         ENDIF
169         CALL tra_zdf_adj_tst( numadt )    ! Vertical mixing
170         IF (lwp) WRITE(numadt,*)
171         CALL tra_nxt_adj_tst( numadt )    ! Tracer fields at next time step
172         IF (lwp) WRITE(numadt,*)
173         CALL istate_init_adj_tst( numadt )
174         IF (lwp) WRITE(numadt,*)
175         CALL tra_sbc_adj_tst( numadt )    ! Surface boundary condition
176         IF (lwp) WRITE(numadt,*)
177         CALL bn2_adj_tst( numadt )        ! Brunt-Vaisala frequency
178         IF (lwp) WRITE(numadt,*)
179         !
180         !-------------- TESTED IN istate_adj_tst ----------------!
181         !CALL zps_hde_adj_tst( numadt )    ! Partial steps: horiz. grad. at bottom level
182         !IF (lwp) WRITE(numadt,*)
183         !--------------------------------------------------------!
184         !
185         ! *** Vertical physics
186         ! --------------------
187         CALL dyn_zdf_adj_tst( numadt )    ! Vertical diffusion
188         IF (lwp) WRITE(numadt,*)
189         CALL dyn_ldf_adj_tst( numadt )    ! Lateral mixing
190         IF (lwp) WRITE(numadt,*)
191         CALL dyn_adv_adj_tst( numadt )    ! Advection (vector or flux form)
192         IF (lwp) WRITE(numadt,*)
193         CALL dyn_hpg_adj_tst( numadt )    ! Horizontal pressure gradient
194         IF (lwp) WRITE(numadt,*)
195         CALL div_cur_adj_tst( numadt )    ! Horizontal divergence and relative vorticity
196         IF (lwp) WRITE(numadt,*)
197         CALL ssh_wzv_adj_tst( numadt )    ! Vertical velocity
198         IF (lwp) WRITE(numadt,*)
199         CALL ssh_nxt_adj_tst( numadt )    ! Sea surface  height time stepping
200         IF (lwp) WRITE(numadt,*)
201         IF ( nn_cla == 1     )  THEN
202            CALL cla_div_adj_tst( numadt )    ! Cross land advection (update hor. divergence)
203            IF (lwp) WRITE(numadt,*)
204         ENDIF
205# if defined key_mpp_mpi
206         CALL lbc_nfd_adj_tst( numadt )
207         IF (lwp) WRITE(numadt,*)
208# endif
209         CALL dyn_vor_adj_tst( numadt )    ! Vorticity term including Coriolis
210         IF (lwp) WRITE(numadt,*)
211         CALL dyn_spg_adj_tst( numadt )    ! Surface pressure gradient
212         IF (lwp) WRITE(numadt,*)
213         IF ( nn_cla == 1 ) THEN
214            CALL cla_dynspg_adj_tst( numadt )    ! Surface pressure gradient
215            IF (lwp) WRITE(numadt,*)
216         END IF
217         CALL dyn_nxt_adj_tst( numadt )    ! Lateral velocity at next time step
218         IF (lwp) WRITE(numadt,*)
219         CALL dyn_bfr_adj_tst( numadt )    ! Surface pressure gradient
220         IF (lwp) WRITE(numadt,*)
221         CALL zdf_bfr_adj_tst( numadt )    ! Surface pressure gradient
222         IF (lwp) WRITE(numadt,*)
223# if defined key_dynspg_flt
224         ! *** Red-Black SOR solver
225         ! ------------
226         CALL sol_sor_adj_tst( numadt )
227         IF (lwp) WRITE(numadt,*)
228# endif
229         CALL flush(numadt)
230         IF (lwp) THEN
231            WRITE(numout,*)
232            WRITE(numout,*) ' tstopt: Finished testing standalone operators'
233            WRITE(numout,*) ' ------'
234            WRITE(numout,*)
235         ENDIF
236
237      ELSEIF (ln_swi_opatam == 1) THEN
238         !
239         ! *** Time-loop operator
240         ! ----------------------
241            CALL stp_adj_tst( numadt )        !Time-stepping
242            IF (lwp) WRITE(numadt,*)
243            CALL flush(numadt)
244      ENDIF
245
246      ! Close output file
247
248      IF (lwp) CLOSE(numadt)
249
250      IF (lwp) THEN
251         WRITE(numout,*)
252         WRITE(numout,*) ' tamtst: Finished testing operators'
253         WRITE(numout,*) ' ------'
254         WRITE(numout,*)
255      ENDIF
256      CALL flush(numout)
257   END SUBROUTINE tam_tst
258   SUBROUTINE tam_tst_init
259      !!----------------------------------------------------------------------
260      !!                     ***  ROUTINE tam_init  ***
261      !!
262      !! ** Purpose :   read tam related namelists and print the variables.
263      !!
264      !! ** input   : - namtst_tam namelist
265      !!              - namtlh namelist
266      !!----------------------------------------------------------------------
267      NAMELIST/namtst_tam/ ln_swi_opatam
268
269
270      ln_swi_opatam  = 0
271
272      REWIND( numnam )              ! Namelist namrun : parameters of the run
273      READ  ( numnam, namtst_tam )
274      IF (lwp) THEN                 ! control print
275         WRITE(numout,*)
276         WRITE(numout,*) 'tam_tst  : Tangent and Adjoint testing'
277         WRITE(numout,*) '~~~~~~~'
278         WRITE(numout,*) '   Namelist namtst_tam'
279         WRITE(numout,*) '      switch for tam testing             ln_swi_opatam  = ', ln_swi_opatam
280      END IF
281   END SUBROUTINE tam_tst_init
282#endif
283END MODULE tamtst
Note: See TracBrowser for help on using the repository browser.