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/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/DIA/diaar5.F90 @ 13519

Last change on this file since 13519 was 13519, checked in by hadcv, 4 years ago

Tiling for DIA routines called by TRA modules

  • Property svn:keywords set to Id
File size: 20.1 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   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   hstr_adv, hstr_ldf
37
38   LOGICAL  :: l_ar5
39     
40   !! * Substitutions
41#  include "do_loop_substitute.h90"
42#  include "domzgr_substitute.h90"
43   !!----------------------------------------------------------------------
44   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
45   !! $Id$
46   !! Software governed by the CeCILL license (see ./LICENSE)
47   !!----------------------------------------------------------------------
48CONTAINS
49
50   FUNCTION dia_ar5_alloc()
51      !!----------------------------------------------------------------------
52      !!                    ***  ROUTINE dia_ar5_alloc  ***
53      !!----------------------------------------------------------------------
54      INTEGER :: dia_ar5_alloc
55      !!----------------------------------------------------------------------
56      !
57      ALLOCATE( thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , &
58         &      hstr_adv(jpi,jpj,jpts,2), hstr_ldf(jpi,jpj,jpts,2), STAT=dia_ar5_alloc )
59      !
60      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
61      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate arrays' )
62      !
63   END FUNCTION dia_ar5_alloc
64
65
66   SUBROUTINE dia_ar5( kt, Kmm )
67      !!----------------------------------------------------------------------
68      !!                    ***  ROUTINE dia_ar5  ***
69      !!
70      !! ** Purpose :   compute and output some AR5 diagnostics
71      !!----------------------------------------------------------------------
72      !
73      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
74      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
75      !
76      INTEGER  ::   ji, jj, jk, iks, ikb                      ! dummy loop arguments
77      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass, zsst
78      REAL(wp) ::   zaw, zbw, zrw
79      !
80      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
81      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z2d, zpe                   ! 2D workspace
82      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: z3d, zrhd, ztpot, zgdept   ! 3D workspace (zgdept: needed to use the substitute)
83      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
84
85      !!--------------------------------------------------------------------
86      IF( ln_timing )   CALL timing_start('dia_ar5')
87 
88      IF( kt == nit000 )     CALL dia_ar5_init
89
90      IF( l_ar5 ) THEN
91         ALLOCATE( zarea_ssh(jpi,jpj), zbotpres(jpi,jpj), z2d(jpi,jpj) )
92         ALLOCATE( zrhd(jpi,jpj,jpk) )
93         ALLOCATE( ztsn(jpi,jpj,jpk,jpts) )
94         zarea_ssh(:,:) = e1e2t(:,:) * ssh(:,:,Kmm)
95      ENDIF
96      !
97      CALL iom_put( 'e2u'      , e2u  (:,:) )
98      CALL iom_put( 'e1v'      , e1v  (:,:) )
99      CALL iom_put( 'areacello', e1e2t(:,:) )
100      !
101      IF( iom_use( 'volcello' ) .OR. iom_use( 'masscello' )  ) THEN 
102         zrhd(:,:,jpk) = 0._wp        ! ocean volume ; rhd is used as workspace
103         DO jk = 1, jpkm1
104            zrhd(:,:,jk) = e1e2t(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
105         END DO
106         DO jk = 1, jpk
107            z3d(:,:,jk) =  rho0 * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
108         END DO
109         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
110         CALL iom_put( 'masscello' , z3d (:,:,:) )   ! ocean mass
111      ENDIF 
112      !
113      IF( iom_use( 'e3tb' ) )  THEN    ! bottom layer thickness
114         DO_2D( 1, 1, 1, 1 )
115            ikb = mbkt(ji,jj)
116            z2d(ji,jj) = e3t(ji,jj,ikb,Kmm)
117         END_2D
118         CALL iom_put( 'e3tb', z2d )
119      ENDIF 
120      !
121      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
122         !                                         ! total volume of liquid seawater
123         zvolssh = glob_sum( 'diaar5', zarea_ssh(:,:) ) 
124         zvol    = vol0 + zvolssh
125     
126         CALL iom_put( 'voltot', zvol               )
127         CALL iom_put( 'sshtot', zvolssh / area_tot )
128         CALL iom_put( 'sshdyn', ssh(:,:,Kmm) - (zvolssh / area_tot) )
129         !
130      ENDIF
131
132      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN   
133         !                     
134         ztsn(:,:,:,jp_tem) = ts(:,:,:,jp_tem,Kmm)                    ! thermosteric ssh
135         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
136         ALLOCATE( zgdept(jpi,jpj,jpk) )
137         DO jk = 1, jpk
138            zgdept(:,:,jk) = gdept(:,:,jk,Kmm)
139         END DO
140         CALL eos( ztsn, zrhd, zgdept)                       ! now in situ density using initial salinity
141         !
142         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
143         DO jk = 1, jpkm1
144            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * zrhd(:,:,jk)
145         END DO
146         IF( ln_linssh ) THEN
147            IF( ln_isfcav ) THEN
148               DO ji = 1, jpi
149                  DO jj = 1, jpj
150                     iks = mikt(ji,jj)
151                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * zrhd(ji,jj,iks) + riceload(ji,jj)
152                  END DO
153               END DO
154            ELSE
155               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
156            END IF
157!!gm
158!!gm   riceload should be added in both ln_linssh=T or F, no?
159!!gm
160         END IF
161         !                                         
162         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
163         zssh_steric = - zarho / area_tot
164         CALL iom_put( 'sshthster', zssh_steric )
165     
166         !                                         ! steric sea surface height
167         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
168         DO jk = 1, jpkm1
169            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * rhd(:,:,jk)
170         END DO
171         IF( ln_linssh ) THEN
172            IF ( ln_isfcav ) THEN
173               DO ji = 1,jpi
174                  DO jj = 1,jpj
175                     iks = mikt(ji,jj)
176                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * rhd(ji,jj,iks) + riceload(ji,jj)
177                  END DO
178               END DO
179            ELSE
180               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * rhd(:,:,1)
181            END IF
182         END IF
183         !   
184         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
185         zssh_steric = - zarho / area_tot
186         CALL iom_put( 'sshsteric', zssh_steric )
187         !                                         ! ocean bottom pressure
188         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
189         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) )
190         CALL iom_put( 'botpres', zbotpres )
191         !
192         DEALLOCATE( zgdept )
193         !
194      ENDIF
195
196      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
197          !                                         ! Mean density anomalie, temperature and salinity
198          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
199          DO_3D( 1, 1, 1, 1, 1, jpkm1 )
200             zztmp = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm)
201             ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * ts(ji,jj,jk,jp_tem,Kmm)
202             ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * ts(ji,jj,jk,jp_sal,Kmm)
203          END_3D
204
205          IF( ln_linssh ) THEN
206            IF( ln_isfcav ) THEN
207               DO ji = 1, jpi
208                  DO jj = 1, jpj
209                     iks = mikt(ji,jj)
210                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_tem,Kmm) 
211                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_sal,Kmm) 
212                  END DO
213               END DO
214            ELSE
215               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * ts(:,:,1,jp_tem,Kmm) 
216               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * ts(:,:,1,jp_sal,Kmm) 
217            END IF
218         ENDIF
219         !
220         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) )
221         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) )
222         zmass = rho0 * ( zarho + zvol )     
223         !
224         CALL iom_put( 'masstot', zmass )
225         CALL iom_put( 'temptot', ztemp / zvol )
226         CALL iom_put( 'saltot' , zsal  / zvol )
227         !
228      ENDIF     
229
230      IF( ln_teos10 ) THEN        ! ! potential temperature (TEOS-10 case)
231         IF( iom_use( 'toce_pot') .OR. iom_use( 'temptot_pot' ) .OR. iom_use( 'sst_pot' )  &
232                                  .OR. iom_use( 'ssttot' ) .OR.  iom_use( 'tosmint_pot' ) ) THEN
233            !
234            ALLOCATE( ztpot(jpi,jpj,jpk) )
235            ztpot(:,:,jpk) = 0._wp
236            DO jk = 1, jpkm1
237               ztpot(:,:,jk) = eos_pt_from_ct( ts(:,:,jk,jp_tem,Kmm), ts(:,:,jk,jp_sal,Kmm) )
238            END DO
239            !
240            CALL iom_put( 'toce_pot', ztpot(:,:,:) )  ! potential temperature (TEOS-10 case)
241            CALL iom_put( 'sst_pot' , ztpot(:,:,1) )  ! surface temperature
242            !
243            IF( iom_use( 'temptot_pot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
244               z2d(:,:) = 0._wp
245               DO jk = 1, jpkm1
246                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t(:,:,jk,Kmm) * ztpot(:,:,jk)
247               END DO
248               ztemp = glob_sum( 'diaar5', z2d(:,:)  ) 
249               CALL iom_put( 'temptot_pot', ztemp / zvol )
250             ENDIF
251             !
252             IF( iom_use( 'ssttot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
253               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  ) 
254               CALL iom_put( 'ssttot', zsst / area_tot )
255             ENDIF
256             ! Vertical integral of temperature
257             IF( iom_use( 'tosmint_pot') ) THEN
258               z2d(:,:) = 0._wp
259               DO_3D( 1, 1, 1, 1, 1, jpkm1 )
260                  z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk)
261               END_3D
262               CALL iom_put( 'tosmint_pot', z2d ) 
263            ENDIF
264            DEALLOCATE( ztpot )
265        ENDIF
266      ELSE       
267         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
268            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * ts(:,:,1,jp_tem,Kmm) )
269            CALL iom_put('ssttot', zsst / area_tot )
270         ENDIF
271      ENDIF
272
273      IF( iom_use( 'tnpeo' )) THEN   
274        ! Work done against stratification by vertical mixing
275        ! Exclude points where rn2 is negative as convection kicks in here and
276        ! work is not being done against stratification
277         ALLOCATE( zpe(jpi,jpj) )
278         zpe(:,:) = 0._wp
279         IF( ln_zdfddm ) THEN
280            DO_3D( 1, 1, 1, 1, 2, jpk )
281               IF( rn2(ji,jj,jk) > 0._wp ) THEN
282                  zrw = ( gdept(ji,jj,jk,Kmm) - gdepw(ji,jj,jk,Kmm) ) / e3w(ji,jj,jk,Kmm)
283                  !
284                  zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
285                  zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
286                  !
287                  zpe(ji, jj) = zpe(ji,jj)   &
288                     &        -  grav * (  avt(ji,jj,jk) * zaw * (ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) )  &
289                     &                   - avs(ji,jj,jk) * zbw * (ts(ji,jj,jk-1,jp_sal,Kmm) - ts(ji,jj,jk,jp_sal,Kmm) ) )
290               ENDIF
291            END_3D
292          ELSE
293            DO_3D( 1, 1, 1, 1, 1, jpk )
294               zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rho0 * e3w(ji,jj,jk,Kmm)
295            END_3D
296         ENDIF
297          CALL iom_put( 'tnpeo', zpe )
298          DEALLOCATE( zpe )
299      ENDIF
300
301      IF( l_ar5 ) THEN
302        DEALLOCATE( zarea_ssh , zbotpres, z2d )
303        DEALLOCATE( ztsn                 )
304      ENDIF
305      !
306      IF( ln_timing )   CALL timing_stop('dia_ar5')
307      !
308   END SUBROUTINE dia_ar5
309
310   ! TEMP: These changes and lbc_lnk not necessary if using XIOS (subdomain support, will not output haloes)
311   SUBROUTINE dia_ar5_hst( ktra, cptr, puflx, pvflx )
312      !!----------------------------------------------------------------------
313      !!                    ***  ROUTINE dia_ar5_htr ***
314      !!----------------------------------------------------------------------
315      !! Wrapper for heat transport calculations
316      !! Called from all advection and/or diffusion routines
317      !!----------------------------------------------------------------------
318      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
319      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
320      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk)    , INTENT(in)   :: puflx  ! u-flux of advection/diffusion
321      REAL(wp), DIMENSION(ST_2D(nn_hls),jpk)    , INTENT(in)   :: pvflx  ! v-flux of advection/diffusion
322      !
323      INTEGER    ::  ji, jj, jk
324
325      IF( cptr /= 'adv' .AND. cptr /= 'ldf' ) RETURN
326      IF( ktra /= jp_tem .AND. ktra /= jp_sal ) RETURN
327
328      IF( cptr == 'adv' ) THEN
329         DO_2D( 0, 0, 0, 0 )
330            hstr_adv(ji,jj,ktra,1) = puflx(ji,jj,1)
331            hstr_adv(ji,jj,ktra,2) = pvflx(ji,jj,1)
332         END_2D
333         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
334            hstr_adv(ji,jj,ktra,1) = hstr_adv(ji,jj,ktra,1) + puflx(ji,jj,jk)
335            hstr_adv(ji,jj,ktra,2) = hstr_adv(ji,jj,ktra,2) + pvflx(ji,jj,jk)
336         END_3D
337      ELSE IF( cptr == 'ldf' ) THEN
338         DO_2D( 0, 0, 0, 0 )
339            hstr_ldf(ji,jj,ktra,1) = puflx(ji,jj,1)
340            hstr_ldf(ji,jj,ktra,2) = pvflx(ji,jj,1)
341         END_2D
342         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
343            hstr_ldf(ji,jj,ktra,1) = hstr_ldf(ji,jj,ktra,1) + puflx(ji,jj,jk)
344            hstr_ldf(ji,jj,ktra,2) = hstr_ldf(ji,jj,ktra,2) + pvflx(ji,jj,jk)
345         END_3D
346      ENDIF
347
348      IF( ntile == 0 .OR. ntile == nijtile ) THEN
349         IF( cptr == 'adv' ) THEN
350            CALL lbc_lnk( 'diaar5', hstr_adv(:,:,ktra,1), 'U', -1.0_wp )
351            IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,1) )  ! advective heat transport in i-direction
352            IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * hstr_adv(:,:,ktra,1) )  ! advective salt transport in i-direction
353            CALL lbc_lnk( 'diaar5', hstr_adv(:,:,ktra,2), 'V', -1.0_wp )
354            IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,2) )  ! advective heat transport in j-direction
355            IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * hstr_adv(:,:,ktra,2) )  ! advective salt transport in j-direction
356         ENDIF
357         IF( cptr == 'ldf' ) THEN
358            CALL lbc_lnk( 'diaar5', hstr_ldf(:,:,ktra,1), 'U', -1.0_wp )
359            IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,1) ) ! diffusive heat transport in i-direction
360            IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * hstr_ldf(:,:,ktra,1) ) ! diffusive salt transport in i-direction
361            CALL lbc_lnk( 'diaar5', hstr_ldf(:,:,ktra,2), 'V', -1.0_wp )
362            IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,2) ) ! diffusive heat transport in j-direction
363            IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * hstr_ldf(:,:,ktra,2) ) ! diffusive salt transport in j-direction
364         ENDIF
365      ENDIF
366   END SUBROUTINE dia_ar5_hst
367
368
369   SUBROUTINE dia_ar5_init
370      !!----------------------------------------------------------------------
371      !!                  ***  ROUTINE dia_ar5_init  ***
372      !!                   
373      !! ** Purpose :   initialization for AR5 diagnostic computation
374      !!----------------------------------------------------------------------
375      INTEGER  ::   inum
376      INTEGER  ::   ik, idep
377      INTEGER  ::   ji, jj, jk  ! dummy loop indices
378      REAL(wp) ::   zztmp 
379      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
380      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
381      !
382      !!----------------------------------------------------------------------
383      !
384      l_ar5 = .FALSE.
385      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
386         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
387         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' ) .OR. &
388         &  iom_use( 'uadv_heattr' ) .OR. iom_use( 'udiff_heattr' ) .OR. &
389         &  iom_use( 'uadv_salttr' ) .OR. iom_use( 'udiff_salttr' ) .OR. &
390         &  iom_use( 'vadv_heattr' ) .OR. iom_use( 'vdiff_heattr' ) .OR. &
391         &  iom_use( 'vadv_salttr' ) .OR. iom_use( 'vdiff_salttr' ) .OR. &
392         &  iom_use( 'rhop' )  ) L_ar5 = .TRUE.
393 
394      IF( l_ar5 ) THEN
395         !
396         !                                      ! allocate dia_ar5 arrays
397         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
398
399         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
400
401         ALLOCATE( zvol0(jpi,jpj) )
402         zvol0 (:,:) = 0._wp
403         thick0(:,:) = 0._wp
404         DO_3D( 1, 1, 1, 1, 1, jpkm1 )
405            idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
406            zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
407            thick0(ji,jj) = thick0(ji,jj) +  idep   
408         END_3D
409         vol0 = glob_sum( 'diaar5', zvol0 )
410         DEALLOCATE( zvol0 )
411
412         IF( iom_use( 'sshthster' ) ) THEN
413            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
414            CALL iom_open ( 'sali_ref_clim_monthly', inum )
415            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,1), 1  )
416            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,2), 12 )
417            CALL iom_close( inum )
418
419            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
420            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
421            IF( ln_zps ) THEN               ! z-coord. partial steps
422               DO_2D( 1, 1, 1, 1 )
423                  ik = mbkt(ji,jj)
424                  IF( ik > 1 ) THEN
425                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
426                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
427                  ENDIF
428               END_2D
429            ENDIF
430            !
431            DEALLOCATE( zsaldta )
432         ENDIF
433         !
434      ENDIF
435      !
436   END SUBROUTINE dia_ar5_init
437
438   !!======================================================================
439END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.