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/2019/dev_r12072_TOP-01_ENHANCE-11_CEthe/src/OCE/DIA – NEMO

source: NEMO/branches/2019/dev_r12072_TOP-01_ENHANCE-11_CEthe/src/OCE/DIA/diaar5.F90 @ 12109

Last change on this file since 12109 was 12109, checked in by cetlod, 4 years ago

check out & merge dev_r11643_ENHANCE-11_CEthe_Shaconemo_diags branch onto dev_r12072_TOP-01_ENHANCE-11_CEthe

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