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/trunk/src/OCE/DIA – NEMO

source: NEMO/trunk/src/OCE/DIA/diaar5.F90 @ 14834

Last change on this file since 14834 was 14834, checked in by hadcv, 3 years ago

#2600: Merge in dev_r14273_HPC-02_Daley_Tiling

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