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

Last change on this file since 12630 was 12630, checked in by jchanut, 4 years ago

Solve #2400 and replace at the same time useless area array by e1e2t

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