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.
diaar5.F90 in NEMO/branches/UKMO/NEMO_4.0.4_CO9_shelf_climate/src/OCE/DIA – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_CO9_shelf_climate/src/OCE/DIA/diaar5.F90 @ 15547

Last change on this file since 15547 was 15333, checked in by hadjt, 3 years ago

Adding Regional Means, but without XIOS or MLD.

Search on
!JT MLD
!JT IOM

in IOM/iom.F90 and DIA/diaregmean.F90

to see the code commented out.

File size: 25.6 KB
Line 
1MODULE diaar5
2   !!======================================================================
3   !!                       ***  MODULE  diaar5  ***
4   !! AR5 diagnostics
5   !!======================================================================
6   !! History :  3.2  !  2009-11  (S. Masson)  Original code
7   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA
8   !!----------------------------------------------------------------------
9   !!   dia_ar5       : AR5 diagnostics
10   !!   dia_ar5_init  : initialisation of AR5 diagnostics
11   !!----------------------------------------------------------------------
12   USE oce            ! ocean dynamics and active tracers
13   USE dom_oce        ! ocean space and time domain
14   USE eosbn2         ! equation of state                (eos_bn2 routine)
15   USE phycst         ! physical constant
16   USE in_out_manager  ! I/O manager
17   USE zdfddm
18   USE zdf_oce
19   !
20   USE lib_mpp        ! distribued memory computing library
21   USE iom            ! I/O manager library
22   USE fldread        ! type FLD_N
23   USE timing         ! preformance summary
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   dia_ar5        ! routine called in step.F90 module
29   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
30   PUBLIC   dia_ar5_hst    ! heat/salt transport
31
32   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
33   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
34   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   thick0       ! ocean thickness (interior domain)
35   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity
36   !JT
37   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   tn0          ! initial temperature
38   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshthster_mat         ! ssh_thermosteric height
39   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshhlster_mat         ! ssh_halosteric height
40   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshsteric_mat         ! ssh_steric height
41   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zbotpres_mat          ! bottom pressure
42
43   !JT
44
45   LOGICAL  :: l_ar5
46     
47   !! * Substitutions
48#  include "vectopt_loop_substitute.h90"
49   !!----------------------------------------------------------------------
50   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
51   !! $Id$
52   !! Software governed by the CeCILL license (see ./LICENSE)
53   !!----------------------------------------------------------------------
54CONTAINS
55
56   FUNCTION dia_ar5_alloc()
57      !!----------------------------------------------------------------------
58      !!                    ***  ROUTINE dia_ar5_alloc  ***
59      !!----------------------------------------------------------------------
60      INTEGER :: dia_ar5_alloc
61      !!----------------------------------------------------------------------
62      !
63      ALLOCATE( thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
64      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
65      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate arrays' )
66
67      !JT
68      ALLOCATE( tn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
69      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
70      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate Temp arrays' )
71
72      ALLOCATE( sshthster_mat(jpi,jpj),sshhlster_mat(jpi,jpj),sshsteric_mat(jpi,jpj), &
73          & zbotpres_mat(jpi,jpj),STAT=dia_ar5_alloc )
74      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
75      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate Temp arrays' )
76
77      !JT
78      !
79   END FUNCTION dia_ar5_alloc
80
81
82   SUBROUTINE dia_ar5( kt )
83      !!----------------------------------------------------------------------
84      !!                    ***  ROUTINE dia_ar5  ***
85      !!
86      !! ** Purpose :   compute and output some AR5 diagnostics
87      !!----------------------------------------------------------------------
88      !
89      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
90      !
91      INTEGER  ::   ji, jj, jk, iks, ikb                      ! dummy loop arguments
92      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass, zsst
93      REAL(wp) ::   zaw, zbw, zrw
94      !
95      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
96      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zpe, z2d                   ! 2D workspace
97      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: zrhd , ztpot               ! 3D workspace
98      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
99
100      !!--------------------------------------------------------------------
101      IF( ln_timing )   CALL timing_start('dia_ar5')
102 
103      IF( kt == nit000 )     CALL dia_ar5_init
104
105      IF( l_ar5 ) THEN
106         ALLOCATE( zarea_ssh(jpi,jpj), zbotpres(jpi,jpj), z2d(jpi,jpj) )
107         ALLOCATE( zrhd(jpi,jpj,jpk) )
108         ALLOCATE( ztsn(jpi,jpj,jpk,jpts) )
109         zarea_ssh(:,:) = e1e2t(:,:) * sshn(:,:)
110      ENDIF
111      !
112      CALL iom_put( 'e2u'      , e2u  (:,:) )
113      CALL iom_put( 'e1v'      , e1v  (:,:) )
114      CALL iom_put( 'areacello', e1e2t(:,:) )
115      !
116      IF( iom_use( 'volcello' ) .OR. iom_use( 'masscello' )  ) THEN 
117         zrhd(:,:,jpk) = 0._wp        ! ocean volume ; rhd is used as workspace
118         DO jk = 1, jpkm1
119            zrhd(:,:,jk) = e1e2t(:,:) * e3t_n(:,:,jk) * tmask(:,:,jk)
120         END DO
121         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
122         CALL iom_put( 'masscello' , rau0 * e3t_n(:,:,:) * tmask(:,:,:) )  ! ocean mass
123      ENDIF 
124      !
125      IF( iom_use( 'e3tb' ) )  THEN    ! bottom layer thickness
126         DO jj = 1, jpj
127            DO ji = 1, jpi
128               ikb = mbkt(ji,jj)
129               z2d(ji,jj) = e3t_n(ji,jj,ikb)
130            END DO
131         END DO
132         CALL iom_put( 'e3tb', z2d )
133      ENDIF 
134      !
135      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
136         !                                         ! total volume of liquid seawater
137         zvolssh = glob_sum( 'diaar5', zarea_ssh(:,:) ) 
138         zvol    = vol0 + zvolssh
139     
140         CALL iom_put( 'voltot', zvol               )
141         CALL iom_put( 'sshtot', zvolssh / area_tot )
142         CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
143         !
144      ENDIF
145
146      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshhlster' )  .OR. iom_use( 'sshsteric' ) .OR. &
147         &  iom_use( 'sshthster_mat' )  .OR. iom_use( 'sshhlster_mat' )  .OR. iom_use( 'sshsteric_mat' )  ) THEN   
148         !                     
149         ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
150         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
151         CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
152         !
153         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
154         DO jk = 1, jpkm1
155            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
156         END DO
157         IF( ln_linssh ) THEN
158            IF( ln_isfcav ) THEN
159               DO ji = 1, jpi
160                  DO jj = 1, jpj
161                     iks = mikt(ji,jj)
162                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
163                  END DO
164               END DO
165            ELSE
166               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
167            END IF
168!!gm
169!!gm   riceload should be added in both ln_linssh=T or F, no?
170!!gm
171         END IF
172         !                                         
173         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
174         zssh_steric = - zarho / area_tot
175         CALL iom_put( 'sshthster', zssh_steric )
176         CALL iom_put( 'sshthster_mat', -zbotpres(:,:) )
177
178
179
180         !JT
181
182 
183         !                     
184         ztsn(:,:,:,jp_tem) = tn0(:,:,:)                    ! halosteric ssh
185         ztsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)                   
186         CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
187         !
188         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
189         DO jk = 1, jpkm1
190            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
191         END DO
192         IF( ln_linssh ) THEN
193            IF( ln_isfcav ) THEN
194               DO ji = 1, jpi
195                  DO jj = 1, jpj
196                     iks = mikt(ji,jj)
197                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
198                  END DO
199               END DO
200            ELSE
201               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
202            END IF
203!!gm
204!!gm   riceload should be added in both ln_linssh=T or F, no?
205!!gm
206         END IF
207         !                                         
208         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
209         zssh_steric = - zarho / area_tot
210         CALL iom_put( 'sshhlster', zssh_steric )
211         CALL iom_put( 'sshhlster_mat', -zbotpres(:,:) )
212
213         !JT
214
215
216
217
218
219
220
221
222
223
224
225     
226         !                                         ! steric sea surface height
227         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
228         DO jk = 1, jpkm1
229            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * rhd(:,:,jk)
230         END DO
231         IF( ln_linssh ) THEN
232            IF ( ln_isfcav ) THEN
233               DO ji = 1,jpi
234                  DO jj = 1,jpj
235                     iks = mikt(ji,jj)
236                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * rhd(ji,jj,iks) + riceload(ji,jj)
237                  END DO
238               END DO
239            ELSE
240               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * rhd(:,:,1)
241            END IF
242         END IF
243         !   
244         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
245         zssh_steric = - zarho / area_tot
246         CALL iom_put( 'sshsteric', zssh_steric )
247         CALL iom_put( 'sshsteric_mat', - zbotpres(:,:) )
248         !                                         ! ocean bottom pressure
249         zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
250         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
251         CALL iom_put( 'botpres', zbotpres )
252         !
253      ENDIF
254
255      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
256          !                                         ! Mean density anomalie, temperature and salinity
257          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
258          DO jk = 1, jpkm1
259             DO jj = 1, jpj
260                DO ji = 1, jpi
261                   zztmp = e1e2t(ji,jj) * e3t_n(ji,jj,jk)
262                   ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * tsn(ji,jj,jk,jp_tem)
263                   ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * tsn(ji,jj,jk,jp_sal)
264                ENDDO
265             ENDDO
266          ENDDO
267
268          IF( ln_linssh ) THEN
269            IF( ln_isfcav ) THEN
270               DO ji = 1, jpi
271                  DO jj = 1, jpj
272                     iks = mikt(ji,jj)
273                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_tem) 
274                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_sal) 
275                  END DO
276               END DO
277            ELSE
278               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * tsn(:,:,1,jp_tem) 
279               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * tsn(:,:,1,jp_sal) 
280            END IF
281         ENDIF
282         !
283         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) )
284         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) )
285         zmass = rau0 * ( zarho + zvol )     
286         !
287         CALL iom_put( 'masstot', zmass )
288         CALL iom_put( 'temptot', ztemp / zvol )
289         CALL iom_put( 'saltot' , zsal  / zvol )
290         !
291      ENDIF     
292
293      IF( ln_teos10 ) THEN        ! ! potential temperature (TEOS-10 case)
294         IF( iom_use( 'toce_pot') .OR. iom_use( 'temptot_pot' ) .OR. iom_use( 'sst_pot' )  &
295                                  .OR. iom_use( 'ssttot' ) .OR.  iom_use( 'tosmint_pot' ) ) THEN
296            !
297            ALLOCATE( ztpot(jpi,jpj,jpk) )
298            ztpot(:,:,jpk) = 0._wp
299            DO jk = 1, jpkm1
300               ztpot(:,:,jk) = eos_pt_from_ct( tsn(:,:,jk,jp_tem), tsn(:,:,jk,jp_sal) )
301            END DO
302            !
303            CALL iom_put( 'toce_pot', ztpot(:,:,:) )  ! potential temperature (TEOS-10 case)
304            CALL iom_put( 'sst_pot' , ztpot(:,:,1) )  ! surface temperature
305            !
306            IF( iom_use( 'temptot_pot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
307               z2d(:,:) = 0._wp
308               DO jk = 1, jpkm1
309                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t_n(:,:,jk) * ztpot(:,:,jk)
310               END DO
311               ztemp = glob_sum( 'diaar5', z2d(:,:)  ) 
312               CALL iom_put( 'temptot_pot', ztemp / zvol )
313             ENDIF
314             !
315             IF( iom_use( 'ssttot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
316               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  ) 
317               CALL iom_put( 'ssttot', zsst / area_tot )
318             ENDIF
319             ! Vertical integral of temperature
320             IF( iom_use( 'tosmint_pot') ) THEN
321               z2d(:,:) = 0._wp
322               DO jk = 1, jpkm1
323                  DO jj = 1, jpj
324                     DO ji = 1, jpi   ! vector opt.
325                        z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t_n(ji,jj,jk) *  ztpot(ji,jj,jk)
326                     END DO
327                  END DO
328               END DO
329               CALL iom_put( 'tosmint_pot', z2d ) 
330            ENDIF
331            DEALLOCATE( ztpot )
332        ENDIF
333      ELSE       
334         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
335            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * tsn(:,:,1,jp_tem) )
336            CALL iom_put('ssttot', zsst / area_tot )
337         ENDIF
338      ENDIF
339
340      IF( iom_use( 'tnpeo' )) THEN   
341        ! Work done against stratification by vertical mixing
342        ! Exclude points where rn2 is negative as convection kicks in here and
343        ! work is not being done against stratification
344         ALLOCATE( zpe(jpi,jpj) )
345         zpe(:,:) = 0._wp
346         IF( ln_zdfddm ) THEN
347            DO jk = 2, jpk
348               DO jj = 1, jpj
349                  DO ji = 1, jpi
350                     IF( rn2(ji,jj,jk) > 0._wp ) THEN
351                        zrw = ( gdept_n(ji,jj,jk) - gdepw_n(ji,jj,jk) ) / e3w_n(ji,jj,jk)
352                        !
353                        zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
354                        zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
355                        !
356                        zpe(ji, jj) = zpe(ji,jj)   &
357                           &        -  grav * (  avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
358                           &                   - avs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
359                     ENDIF
360                  END DO
361               END DO
362             END DO
363          ELSE
364            DO jk = 1, jpk
365               DO ji = 1, jpi
366                  DO jj = 1, jpj
367                     zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rau0 * e3w_n(ji,jj,jk)
368                  END DO
369               END DO
370            END DO
371         ENDIF
372          CALL iom_put( 'tnpeo', zpe )
373          DEALLOCATE( zpe )
374      ENDIF
375
376      IF( l_ar5 ) THEN
377        DEALLOCATE( zarea_ssh , zbotpres, z2d )
378        DEALLOCATE( zrhd      )
379        DEALLOCATE( ztsn      )
380      ENDIF
381      !
382      IF( ln_timing )   CALL timing_stop('dia_ar5')
383      !
384   END SUBROUTINE dia_ar5
385
386
387   SUBROUTINE dia_ar5_hst( ktra, cptr, pua, pva ) 
388      !!----------------------------------------------------------------------
389      !!                    ***  ROUTINE dia_ar5_htr ***
390      !!----------------------------------------------------------------------
391      !! Wrapper for heat transport calculations
392      !! Called from all advection and/or diffusion routines
393      !!----------------------------------------------------------------------
394      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
395      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
396      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pua   ! 3D input array of advection/diffusion
397      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pva   ! 3D input array of advection/diffusion
398      !
399      INTEGER    ::  ji, jj, jk
400      REAL(wp), DIMENSION(jpi,jpj)  :: z2d
401
402   
403      z2d(:,:) = pua(:,:,1) 
404      DO jk = 1, jpkm1
405         DO jj = 2, jpjm1
406            DO ji = fs_2, fs_jpim1   ! vector opt.
407               z2d(ji,jj) = z2d(ji,jj) + pua(ji,jj,jk) 
408            END DO
409         END DO
410       END DO
411       CALL lbc_lnk( 'diaar5', z2d, 'U', -1. )
412       IF( cptr == 'adv' ) THEN
413          IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in i-direction
414          IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rau0     * z2d )  ! advective salt transport in i-direction
415       ENDIF
416       IF( cptr == 'ldf' ) THEN
417          IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in i-direction
418          IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rau0     * z2d ) ! diffusive salt transport in i-direction
419       ENDIF
420       !
421       z2d(:,:) = pva(:,:,1) 
422       DO jk = 1, jpkm1
423          DO jj = 2, jpjm1
424             DO ji = fs_2, fs_jpim1   ! vector opt.
425                z2d(ji,jj) = z2d(ji,jj) + pva(ji,jj,jk) 
426             END DO
427          END DO
428       END DO
429       CALL lbc_lnk( 'diaar5', z2d, 'V', -1. )
430       IF( cptr == 'adv' ) THEN
431          IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in j-direction
432          IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rau0     * z2d )  ! advective salt transport in j-direction
433       ENDIF
434       IF( cptr == 'ldf' ) THEN
435          IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in j-direction
436          IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rau0     * z2d ) ! diffusive salt transport in j-direction
437       ENDIF
438         
439   END SUBROUTINE dia_ar5_hst
440
441
442   SUBROUTINE dia_ar5_init
443      !!----------------------------------------------------------------------
444      !!                  ***  ROUTINE dia_ar5_init  ***
445      !!                   
446      !! ** Purpose :   initialization for AR5 diagnostic computation
447      !!----------------------------------------------------------------------
448      INTEGER  ::   inum
449      INTEGER  ::   ik, idep
450      INTEGER  ::   ji, jj, jk  ! dummy loop indices
451      REAL(wp) ::   zztmp 
452      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
453      !JT
454      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   ztemdta   ! Jan/Dec levitus salinity
455      !JT
456      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
457
458
459      !JT
460      !                                  !!* namtsd  namelist : Temperature & Salinity Data *
461      LOGICAL ::   ln_tsd_init   !: T & S data flag
462      LOGICAL ::   ln_tsd_dmp    !: internal damping toward input data flag
463      INTEGER ::   ios, ierr0, ierr1, ierr2, ierr3   ! local integers
464
465      CHARACTER(len=100)            ::   cn_dir          ! Root directory for location of ssr files
466      ! TYPE(FLD_N), DIMENSION( jpts) ::   slf_i           ! array of namelist informations on the fields to read
467      TYPE(FLD_N)                   ::   sn_tem, sn_sal
468      !JT
469      !
470      !!----------------------------------------------------------------------
471      !
472      l_ar5 = .FALSE.
473      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
474         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
475         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' ) .OR. &
476         &  iom_use( 'rhop' )    .OR. iom_use( 'sshhlster' )  .OR. &
477         &  iom_use( 'sshthster_mat' )  .OR. iom_use( 'sshsteric_mat' ) .OR. iom_use( 'sshhlster_mat' )   ) L_ar5 = .TRUE.
478 
479      IF( l_ar5 ) THEN
480         !
481         !                                      ! allocate dia_ar5 arrays
482         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
483
484         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
485
486         ALLOCATE( zvol0(jpi,jpj) )
487         zvol0 (:,:) = 0._wp
488         thick0(:,:) = 0._wp
489         DO jk = 1, jpkm1
490            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
491               DO ji = 1, jpi
492                  idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
493                  zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
494                  thick0(ji,jj) = thick0(ji,jj) +  idep   
495               END DO
496            END DO
497         END DO
498         vol0 = glob_sum( 'diaar5', zvol0 )
499         DEALLOCATE( zvol0 )
500
501         !JT
502         IF( iom_use( 'sshthster' ) .OR. iom_use( 'sshthster_mat' ) .OR. & 
503            & iom_use( 'sshhlster' ) .OR. iom_use( 'sshhlster_mat' )  ) THEN
504
505            NAMELIST/namtsd/   ln_tsd_init, ln_tsd_dmp, cn_dir, sn_tem, sn_sal
506            !!----------------------------------------------------------------------
507            !
508            !  Initialisation
509            ierr0 = 0  ;  ierr1 = 0  ;  ierr2 = 0  ;  ierr3 = 0
510            !
511            REWIND( numnam_ref )              ! Namelist namtsd in reference namelist :
512            READ  ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901)
513901         IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtsd in reference namelist' )
514            REWIND( numnam_cfg )              ! Namelist namtsd in configuration namelist : Parameters of the run
515            READ  ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 )
516902         IF( ios >  0 )   CALL ctl_nam ( ios , 'namtsd in configuration namelist' )
517            IF(lwm) WRITE ( numond, namtsd )
518             
519            IF(lwp) THEN                  ! control print
520                 WRITE(numout,*)
521                 WRITE(numout,*) 'dia_ar5_init : Temperature & Salinity data from namtsd '
522                 WRITE(numout,*) '~~~~~~~~~~~~ '
523                 WRITE(numout,*) '   Namelist namtsd'
524                 WRITE(numout,*) '      T data   ln_tsd_init = ', ln_tsd_init
525                 WRITE(numout,*) '      damping of ocean T & S toward T &S input data        ln_tsd_dmp  = ', ln_tsd_dmp
526                 WRITE(numout,*)
527            ENDIF
528             
529          ENDIF
530          IF( iom_use( 'sshthster' ) .OR. iom_use( 'sshthster_mat' )   ) THEN
531             
532              ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
533              CALL iom_open ( TRIM( cn_dir )//TRIM(sn_sal%clname), inum )
534              CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,1), 1  )
535              CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,2), 12 )
536              CALL iom_close( inum )
537             
538              sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
539              sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
540             
541              IF( ln_zps ) THEN               ! z-coord. partial steps
542                 DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
543                    DO ji = 1, jpi
544                       ik = mbkt(ji,jj)
545                       IF( ik > 1 ) THEN
546                          zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
547                          sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
548                       ENDIF
549                    END DO
550                 END DO
551              ENDIF
552              !
553              DEALLOCATE( zsaldta )
554          ENDIF
555          IF( iom_use( 'sshhlster' ) .OR. iom_use( 'sshhlster_mat' )   ) THEN
556         
557              ALLOCATE( ztemdta(jpi,jpj,jpk,jpts) )
558              CALL iom_open ( TRIM( cn_dir )//TRIM(sn_tem%clname), inum )
559              CALL iom_get  ( inum, jpdom_data, TRIM(sn_tem%clvar), ztemdta(:,:,:,1), 1  )
560              CALL iom_get  ( inum, jpdom_data, TRIM(sn_tem%clvar), ztemdta(:,:,:,2), 12 )
561              CALL iom_close( inum )
562
563              tn0(:,:,:) = 0.5_wp * ( ztemdta(:,:,:,1) + ztemdta(:,:,:,2) )       
564              tn0(:,:,:) = tn0(:,:,:) * tmask(:,:,:)
565
566              IF( ln_zps ) THEN               ! z-coord. partial steps
567                 DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
568                    DO ji = 1, jpi
569                       ik = mbkt(ji,jj)
570                       IF( ik > 1 ) THEN
571                          zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
572                          tn0(ji,jj,ik) = ( 1._wp - zztmp ) * tn0(ji,jj,ik) + zztmp * tn0(ji,jj,ik-1)
573                       ENDIF
574                    END DO
575                 END DO
576              ENDIF
577              !
578              DEALLOCATE( ztemdta )
579          ENDIF
580     
581      ENDIF
582         !JT
583         !
584      !
585   END SUBROUTINE dia_ar5_init
586
587   !!======================================================================
588END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.