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/r4.0-HEAD_r12713_dan_test_clems_branch/src/OCE/DIA – NEMO

source: NEMO/branches/UKMO/r4.0-HEAD_r12713_dan_test_clems_branch/src/OCE/DIA/diaar5.F90

Last change on this file was 12801, checked in by dancopsey, 4 years ago

Clear SVN keywords

File size: 19.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
40#  include "vectopt_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      !
[12631]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
[1756]63   SUBROUTINE dia_ar5( kt )
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
[2715]71      !
[12276]72      INTEGER  ::   ji, jj, jk, iks, ikb                      ! dummy loop arguments
73      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass, zsst
[7646]74      REAL(wp) ::   zaw, zbw, zrw
[3294]75      !
[9125]76      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
[12276]77      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zpe, z2d                   ! 2D workspace
78      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: zrhd , zrhop, ztpot   ! 3D workspace
[9125]79      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
80
[1756]81      !!--------------------------------------------------------------------
[9124]82      IF( ln_timing )   CALL timing_start('dia_ar5')
[3294]83 
[7646]84      IF( kt == nit000 )     CALL dia_ar5_init
[1756]85
[9125]86      IF( l_ar5 ) THEN
[12276]87         ALLOCATE( zarea_ssh(jpi,jpj), zbotpres(jpi,jpj), z2d(jpi,jpj) )
[9125]88         ALLOCATE( zrhd(jpi,jpj,jpk) , zrhop(jpi,jpj,jpk) )
89         ALLOCATE( ztsn(jpi,jpj,jpk,jpts) )
[12631]90         zarea_ssh(:,:) = e1e2t(:,:) * sshn(:,:)
[7646]91      ENDIF
92      !
[12631]93      CALL iom_put( 'e2u'      , e2u  (:,:) )
94      CALL iom_put( 'e1v'      , e1v  (:,:) )
95      CALL iom_put( 'areacello', e1e2t(:,:) )
[12276]96      !
97      IF( iom_use( 'volcello' ) .OR. iom_use( 'masscello' )  ) THEN 
98         zrhd(:,:,jpk) = 0._wp        ! ocean volume ; rhd is used as workspace
99         DO jk = 1, jpkm1
[12631]100            zrhd(:,:,jk) = e1e2t(:,:) * e3t_n(:,:,jk) * tmask(:,:,jk)
[12276]101         END DO
102         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
103         CALL iom_put( 'masscello' , rau0 * e3t_n(:,:,:) * tmask(:,:,:) )  ! ocean mass
104      ENDIF 
105      !
106      IF( iom_use( 'e3tb' ) )  THEN    ! bottom layer thickness
107         DO jj = 1, jpj
108            DO ji = 1, jpi
109               ikb = mbkt(ji,jj)
110               z2d(ji,jj) = e3t_n(ji,jj,ikb)
111            END DO
112         END DO
113         CALL iom_put( 'e3tb', z2d )
114      ENDIF 
115      !
[7646]116      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
117         !                                         ! total volume of liquid seawater
[12276]118         zvolssh = glob_sum( 'diaar5', zarea_ssh(:,:) ) 
119         zvol    = vol0 + zvolssh
[1756]120     
[7646]121         CALL iom_put( 'voltot', zvol               )
122         CALL iom_put( 'sshtot', zvolssh / area_tot )
123         CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
124         !
125      ENDIF
[1756]126
[7646]127      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN   
128         !                     
[7753]129         ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
130         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
[7646]131         CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
132         !
[7753]133         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]134         DO jk = 1, jpkm1
[7753]135            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
[7646]136         END DO
137         IF( ln_linssh ) THEN
138            IF( ln_isfcav ) THEN
139               DO ji = 1, jpi
140                  DO jj = 1, jpj
[12276]141                     iks = mikt(ji,jj)
142                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
[7646]143                  END DO
[5120]144               END DO
[7646]145            ELSE
[7753]146               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
[7646]147            END IF
[6140]148!!gm
149!!gm   riceload should be added in both ln_linssh=T or F, no?
150!!gm
[7646]151         END IF
152         !                                         
[12631]153         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
[7646]154         zssh_steric = - zarho / area_tot
155         CALL iom_put( 'sshthster', zssh_steric )
[1756]156     
[7646]157         !                                         ! steric sea surface height
158         CALL eos( tsn, zrhd, zrhop, gdept_n(:,:,:) )                 ! now in situ and potential density
[7753]159         zrhop(:,:,jpk) = 0._wp
[7646]160         CALL iom_put( 'rhop', zrhop )
161         !
[7753]162         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]163         DO jk = 1, jpkm1
[7753]164            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
[7646]165         END DO
166         IF( ln_linssh ) THEN
167            IF ( ln_isfcav ) THEN
168               DO ji = 1,jpi
169                  DO jj = 1,jpj
[12276]170                     iks = mikt(ji,jj)
171                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
[7646]172                  END DO
[5120]173               END DO
[7646]174            ELSE
[7753]175               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
[7646]176            END IF
[5120]177         END IF
[7646]178         !   
[12631]179         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
[7646]180         zssh_steric = - zarho / area_tot
181         CALL iom_put( 'sshsteric', zssh_steric )
182         !                                         ! ocean bottom pressure
183         zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
[7753]184         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
[7646]185         CALL iom_put( 'botpres', zbotpres )
186         !
187      ENDIF
[1756]188
[7646]189      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
[12276]190          !                                         ! Mean density anomalie, temperature and salinity
191          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
192          DO jk = 1, jpkm1
193             DO jj = 1, jpj
194                DO ji = 1, jpi
[12631]195                   zztmp = e1e2t(ji,jj) * e3t_n(ji,jj,jk)
[12276]196                   ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * tsn(ji,jj,jk,jp_tem)
197                   ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * tsn(ji,jj,jk,jp_sal)
198                ENDDO
199             ENDDO
200          ENDDO
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)
207                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_tem) 
208                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_sal) 
[7646]209                  END DO
[5120]210               END DO
[7646]211            ELSE
[12276]212               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * tsn(:,:,1,jp_tem) 
213               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * tsn(:,:,1,jp_sal) 
[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
234               ztpot(:,:,jk) = eos_pt_from_ct( tsn(:,:,jk,jp_tem), tsn(:,:,jk,jp_sal) )
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
[12631]243                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t_n(:,:,jk) * 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
[12631]250               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  ) 
[12276]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
256               DO jk = 1, jpkm1
257                  DO jj = 1, jpj
258                     DO ji = 1, jpi   ! vector opt.
259                        z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t_n(ji,jj,jk) *  ztpot(ji,jj,jk)
260                     END DO
261                  END DO
262               END DO
263               CALL iom_put( 'tosmint_pot', z2d ) 
264            ENDIF
265            DEALLOCATE( ztpot )
266        ENDIF
267      ELSE       
268         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
[12631]269            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * tsn(:,:,1,jp_tem) )
[12276]270            CALL iom_put('ssttot', zsst / area_tot )
271         ENDIF
[1756]272      ENDIF
[7646]273
274      IF( iom_use( 'tnpeo' )) THEN   
[12276]275        ! Work done against stratification by vertical mixing
276        ! Exclude points where rn2 is negative as convection kicks in here and
277        ! work is not being done against stratification
[9125]278         ALLOCATE( zpe(jpi,jpj) )
[8078]279         zpe(:,:) = 0._wp
[9019]280         IF( ln_zdfddm ) THEN
[8078]281            DO jk = 2, jpk
282               DO jj = 1, jpj
283                  DO ji = 1, jpi
284                     IF( rn2(ji,jj,jk) > 0._wp ) THEN
[12276]285                        zrw = ( gdept_n(ji,jj,jk) - gdepw_n(ji,jj,jk) ) / e3w_n(ji,jj,jk)
[8078]286                        !
287                        zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
288                        zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
289                        !
[12276]290                        zpe(ji, jj) = zpe(ji,jj)   &
[9019]291                           &        -  grav * (  avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
292                           &                   - avs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
[8078]293                     ENDIF
294                  END DO
295               END DO
296             END DO
[7646]297          ELSE
[8078]298            DO jk = 1, jpk
299               DO ji = 1, jpi
300                  DO jj = 1, jpj
[12276]301                     zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rau0 * e3w_n(ji,jj,jk)
[8078]302                  END DO
303               END DO
304            END DO
305         ENDIF
[9019]306          CALL iom_put( 'tnpeo', zpe )
[9125]307          DEALLOCATE( zpe )
[7646]308      ENDIF
[9125]309
[7646]310      IF( l_ar5 ) THEN
[12276]311        DEALLOCATE( zarea_ssh , zbotpres, z2d )
[9125]312        DEALLOCATE( zrhd      , zrhop    )
313        DEALLOCATE( ztsn                 )
[7646]314      ENDIF
[1756]315      !
[9124]316      IF( ln_timing )   CALL timing_stop('dia_ar5')
[3294]317      !
[1756]318   END SUBROUTINE dia_ar5
319
[9124]320
[7646]321   SUBROUTINE dia_ar5_hst( ktra, cptr, pua, pva ) 
322      !!----------------------------------------------------------------------
323      !!                    ***  ROUTINE dia_ar5_htr ***
324      !!----------------------------------------------------------------------
325      !! Wrapper for heat transport calculations
326      !! Called from all advection and/or diffusion routines
327      !!----------------------------------------------------------------------
328      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
329      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
330      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pua   ! 3D input array of advection/diffusion
331      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pva   ! 3D input array of advection/diffusion
332      !
333      INTEGER    ::  ji, jj, jk
[9125]334      REAL(wp), DIMENSION(jpi,jpj)  :: z2d
[1756]335
[7646]336   
337      z2d(:,:) = pua(:,:,1) 
338      DO jk = 1, jpkm1
339         DO jj = 2, jpjm1
340            DO ji = fs_2, fs_jpim1   ! vector opt.
341               z2d(ji,jj) = z2d(ji,jj) + pua(ji,jj,jk) 
342            END DO
343         END DO
344       END DO
[10425]345       CALL lbc_lnk( 'diaar5', z2d, 'U', -1. )
[7646]346       IF( cptr == 'adv' ) THEN
[12276]347          IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in i-direction
348          IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rau0     * z2d )  ! advective salt transport in i-direction
[7646]349       ENDIF
350       IF( cptr == 'ldf' ) THEN
[12276]351          IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in i-direction
352          IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rau0     * z2d ) ! diffusive salt transport in i-direction
[7646]353       ENDIF
354       !
355       z2d(:,:) = pva(:,:,1) 
356       DO jk = 1, jpkm1
357          DO jj = 2, jpjm1
358             DO ji = fs_2, fs_jpim1   ! vector opt.
359                z2d(ji,jj) = z2d(ji,jj) + pva(ji,jj,jk) 
360             END DO
361          END DO
362       END DO
[10425]363       CALL lbc_lnk( 'diaar5', z2d, 'V', -1. )
[7646]364       IF( cptr == 'adv' ) THEN
[12276]365          IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in j-direction
366          IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rau0     * z2d )  ! advective salt transport in j-direction
[7646]367       ENDIF
368       IF( cptr == 'ldf' ) THEN
[12276]369          IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in j-direction
370          IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rau0     * z2d ) ! diffusive salt transport in j-direction
[7646]371       ENDIF
372         
373   END SUBROUTINE dia_ar5_hst
374
375
[1756]376   SUBROUTINE dia_ar5_init
377      !!----------------------------------------------------------------------
378      !!                  ***  ROUTINE dia_ar5_init  ***
379      !!                   
[2528]380      !! ** Purpose :   initialization for AR5 diagnostic computation
[1756]381      !!----------------------------------------------------------------------
382      INTEGER  ::   inum
[12276]383      INTEGER  ::   ik, idep
[1756]384      INTEGER  ::   ji, jj, jk  ! dummy loop indices
[7753]385      REAL(wp) ::   zztmp 
[9125]386      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
[12276]387      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
[5253]388      !
[1756]389      !!----------------------------------------------------------------------
390      !
[7646]391      l_ar5 = .FALSE.
392      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
393         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
394         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) L_ar5 = .TRUE.
395 
396      IF( l_ar5 ) THEN
397         !
398         !                                      ! allocate dia_ar5 arrays
399         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
[2715]400
[12631]401         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
[1756]402
[12276]403         ALLOCATE( zvol0(jpi,jpj) )
404         zvol0 (:,:) = 0._wp
[7753]405         thick0(:,:) = 0._wp
[7646]406         DO jk = 1, jpkm1
[12276]407            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
408               DO ji = 1, jpi
409                  idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
[12631]410                  zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
[12276]411                  thick0(ji,jj) = thick0(ji,jj) +  idep   
412               END DO
413            END DO
[7646]414         END DO
[12276]415         vol0 = glob_sum( 'diaar5', zvol0 )
416         DEALLOCATE( zvol0 )
[5253]417
[9125]418         IF( iom_use( 'sshthster' ) ) THEN
[12276]419            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
[9125]420            CALL iom_open ( 'sali_ref_clim_monthly', inum )
421            CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
422            CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
423            CALL iom_close( inum )
[6665]424
[9125]425            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
426            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
427            IF( ln_zps ) THEN               ! z-coord. partial steps
428               DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
429                  DO ji = 1, jpi
430                     ik = mbkt(ji,jj)
431                     IF( ik > 1 ) THEN
432                        zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
433                        sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
434                     ENDIF
435                  END DO
[7646]436               END DO
[9125]437            ENDIF
438            !
439            DEALLOCATE( zsaldta )
[7646]440         ENDIF
441         !
[1756]442      ENDIF
443      !
444   END SUBROUTINE dia_ar5_init
445
446   !!======================================================================
447END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.