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_r12377_KERNEL-06_techene_e3/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/DIA/diaar5.F90 @ 12680

Last change on this file since 12680 was 12680, checked in by techene, 4 years ago

dynatfQCO.F90, stepLF.F90 : fixed (remove pe3. from dyn_atf_qco input arguments), all : remove e3. tables and include gurvan's feedbacks

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