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_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE/DIA/diaar5.F90 @ 14037

Last change on this file since 14037 was 14037, checked in by ayoung, 3 years ago

Updated to trunk at 14020. Sette tests passed with change of results for configurations with non-linear ssh. Ticket #2506.

  • Property svn:keywords set to Id
File size: 19.9 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
[14037]36   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   hstr_adv, hstr_ldf
[7646]37
38   LOGICAL  :: l_ar5
[1756]39     
[7646]40   !! * Substitutions
[12377]41#  include "do_loop_substitute.h90"
[13237]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      !
[14037]57      ALLOCATE( thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , &
58         &      hstr_adv(jpi,jpj,jpts,2), hstr_ldf(jpi,jpj,jpts,2), STAT=dia_ar5_alloc )
[2715]59      !
[10425]60      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
61      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate arrays' )
[2715]62      !
63   END FUNCTION dia_ar5_alloc
64
65
[12377]66   SUBROUTINE dia_ar5( kt, Kmm )
[1756]67      !!----------------------------------------------------------------------
68      !!                    ***  ROUTINE dia_ar5  ***
69      !!
[2528]70      !! ** Purpose :   compute and output some AR5 diagnostics
[1756]71      !!----------------------------------------------------------------------
[2715]72      !
[1756]73      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
[12377]74      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
[2715]75      !
[12276]76      INTEGER  ::   ji, jj, jk, iks, ikb                      ! dummy loop arguments
77      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass, zsst
[7646]78      REAL(wp) ::   zaw, zbw, zrw
[3294]79      !
[9125]80      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
[13237]81      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z2d, zpe                   ! 2D workspace
82      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: z3d, zrhd, ztpot, zgdept   ! 3D workspace (zgdept: needed to use the substitute)
[9125]83      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
84
[1756]85      !!--------------------------------------------------------------------
[9124]86      IF( ln_timing )   CALL timing_start('dia_ar5')
[3294]87 
[7646]88      IF( kt == nit000 )     CALL dia_ar5_init
[1756]89
[9125]90      IF( l_ar5 ) THEN
[12276]91         ALLOCATE( zarea_ssh(jpi,jpj), zbotpres(jpi,jpj), z2d(jpi,jpj) )
[13089]92         ALLOCATE( zrhd(jpi,jpj,jpk) )
[9125]93         ALLOCATE( ztsn(jpi,jpj,jpk,jpts) )
[12630]94         zarea_ssh(:,:) = e1e2t(:,:) * ssh(:,:,Kmm)
[7646]95      ENDIF
96      !
[12630]97      CALL iom_put( 'e2u'      , e2u  (:,:) )
98      CALL iom_put( 'e1v'      , e1v  (:,:) )
99      CALL iom_put( 'areacello', e1e2t(:,:) )
[12276]100      !
101      IF( iom_use( 'volcello' ) .OR. iom_use( 'masscello' )  ) THEN 
102         zrhd(:,:,jpk) = 0._wp        ! ocean volume ; rhd is used as workspace
103         DO jk = 1, jpkm1
[12630]104            zrhd(:,:,jk) = e1e2t(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
[12276]105         END DO
[13237]106         DO jk = 1, jpk
107            z3d(:,:,jk) =  rho0 * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
108         END DO
[12276]109         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
[13237]110         CALL iom_put( 'masscello' , z3d (:,:,:) )   ! ocean mass
[12276]111      ENDIF 
112      !
113      IF( iom_use( 'e3tb' ) )  THEN    ! bottom layer thickness
[13295]114         DO_2D( 1, 1, 1, 1 )
[12377]115            ikb = mbkt(ji,jj)
116            z2d(ji,jj) = e3t(ji,jj,ikb,Kmm)
117         END_2D
[12276]118         CALL iom_put( 'e3tb', z2d )
119      ENDIF 
120      !
[7646]121      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
122         !                                         ! total volume of liquid seawater
[12276]123         zvolssh = glob_sum( 'diaar5', zarea_ssh(:,:) ) 
124         zvol    = vol0 + zvolssh
[1756]125     
[7646]126         CALL iom_put( 'voltot', zvol               )
127         CALL iom_put( 'sshtot', zvolssh / area_tot )
[12377]128         CALL iom_put( 'sshdyn', ssh(:,:,Kmm) - (zvolssh / area_tot) )
[7646]129         !
130      ENDIF
[1756]131
[7646]132      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN   
133         !                     
[12377]134         ztsn(:,:,:,jp_tem) = ts(:,:,:,jp_tem,Kmm)                    ! thermosteric ssh
[7753]135         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
[13237]136         ALLOCATE( zgdept(jpi,jpj,jpk) )
137         DO jk = 1, jpk
138            zgdept(:,:,jk) = gdept(:,:,jk,Kmm)
139         END DO
140         CALL eos( ztsn, zrhd, zgdept)                       ! now in situ density using initial salinity
[7646]141         !
[7753]142         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]143         DO jk = 1, jpkm1
[12377]144            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * zrhd(:,:,jk)
[7646]145         END DO
146         IF( ln_linssh ) THEN
147            IF( ln_isfcav ) THEN
[14037]148               DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
149                  iks = mikt(ji,jj)
150                  zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * zrhd(ji,jj,iks) + riceload(ji,jj)
151               END_2D
[7646]152            ELSE
[12377]153               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
[7646]154            END IF
[6140]155!!gm
156!!gm   riceload should be added in both ln_linssh=T or F, no?
157!!gm
[7646]158         END IF
159         !                                         
[12630]160         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
[7646]161         zssh_steric = - zarho / area_tot
162         CALL iom_put( 'sshthster', zssh_steric )
[1756]163     
[7646]164         !                                         ! steric sea surface height
[7753]165         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]166         DO jk = 1, jpkm1
[13089]167            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * rhd(:,:,jk)
[7646]168         END DO
169         IF( ln_linssh ) THEN
170            IF ( ln_isfcav ) THEN
171               DO ji = 1,jpi
172                  DO jj = 1,jpj
[12276]173                     iks = mikt(ji,jj)
[13089]174                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * rhd(ji,jj,iks) + riceload(ji,jj)
[7646]175                  END DO
[5120]176               END DO
[7646]177            ELSE
[13089]178               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * rhd(:,:,1)
[7646]179            END IF
[5120]180         END IF
[7646]181         !   
[12630]182         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
[7646]183         zssh_steric = - zarho / area_tot
184         CALL iom_put( 'sshsteric', zssh_steric )
185         !                                         ! ocean bottom pressure
[12489]186         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
[12377]187         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) )
[7646]188         CALL iom_put( 'botpres', zbotpres )
189         !
[13237]190         DEALLOCATE( zgdept )
191         !
[7646]192      ENDIF
[1756]193
[7646]194      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
[12276]195          !                                         ! Mean density anomalie, temperature and salinity
196          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
[13295]197          DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[12630]198             zztmp = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm)
[12377]199             ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * ts(ji,jj,jk,jp_tem,Kmm)
200             ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * ts(ji,jj,jk,jp_sal,Kmm)
201          END_3D
[12276]202
203          IF( ln_linssh ) THEN
[7646]204            IF( ln_isfcav ) THEN
205               DO ji = 1, jpi
206                  DO jj = 1, jpj
[12276]207                     iks = mikt(ji,jj)
[12377]208                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_tem,Kmm) 
209                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_sal,Kmm) 
[7646]210                  END DO
[5120]211               END DO
[7646]212            ELSE
[12377]213               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * ts(:,:,1,jp_tem,Kmm) 
214               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * ts(:,:,1,jp_sal,Kmm) 
[7646]215            END IF
216         ENDIF
217         !
[12276]218         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) )
219         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) )
[12489]220         zmass = rho0 * ( zarho + zvol )     
[7646]221         !
222         CALL iom_put( 'masstot', zmass )
[12276]223         CALL iom_put( 'temptot', ztemp / zvol )
224         CALL iom_put( 'saltot' , zsal  / zvol )
[7646]225         !
[12276]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
[12377]235               ztpot(:,:,jk) = eos_pt_from_ct( ts(:,:,jk,jp_tem,Kmm), ts(:,:,jk,jp_sal,Kmm) )
[12276]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
[12630]244                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t(:,:,jk,Kmm) * ztpot(:,:,jk)
[12276]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
[12630]251               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  ) 
[12276]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
[13295]257               DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[12489]258                  z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk)
[12377]259               END_3D
[12276]260               CALL iom_put( 'tosmint_pot', z2d ) 
261            ENDIF
262            DEALLOCATE( ztpot )
263        ENDIF
264      ELSE       
265         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
[12630]266            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * ts(:,:,1,jp_tem,Kmm) )
[12276]267            CALL iom_put('ssttot', zsst / area_tot )
268         ENDIF
[1756]269      ENDIF
[7646]270
271      IF( iom_use( 'tnpeo' )) THEN   
[12276]272        ! Work done against stratification by vertical mixing
273        ! Exclude points where rn2 is negative as convection kicks in here and
274        ! work is not being done against stratification
[9125]275         ALLOCATE( zpe(jpi,jpj) )
[8078]276         zpe(:,:) = 0._wp
[9019]277         IF( ln_zdfddm ) THEN
[13295]278            DO_3D( 1, 1, 1, 1, 2, jpk )
[12377]279               IF( rn2(ji,jj,jk) > 0._wp ) THEN
280                  zrw = ( gdept(ji,jj,jk,Kmm) - gdepw(ji,jj,jk,Kmm) ) / e3w(ji,jj,jk,Kmm)
281                  !
282                  zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
283                  zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
284                  !
285                  zpe(ji, jj) = zpe(ji,jj)   &
286                     &        -  grav * (  avt(ji,jj,jk) * zaw * (ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) )  &
287                     &                   - avs(ji,jj,jk) * zbw * (ts(ji,jj,jk-1,jp_sal,Kmm) - ts(ji,jj,jk,jp_sal,Kmm) ) )
288               ENDIF
289            END_3D
[7646]290          ELSE
[13295]291            DO_3D( 1, 1, 1, 1, 1, jpk )
[12489]292               zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rho0 * e3w(ji,jj,jk,Kmm)
[12377]293            END_3D
[8078]294         ENDIF
[9019]295          CALL iom_put( 'tnpeo', zpe )
[9125]296          DEALLOCATE( zpe )
[7646]297      ENDIF
[9125]298
[7646]299      IF( l_ar5 ) THEN
[12276]300        DEALLOCATE( zarea_ssh , zbotpres, z2d )
[9125]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
[14037]308   ! TEMP: [tiling] These changes not necessary if using XIOS (subdomain support, will not output haloes)
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'
[14037]318      REAL(wp), DIMENSION(A2D(nn_hls),jpk)    , INTENT(in)   :: puflx  ! u-flux of advection/diffusion
319      REAL(wp), DIMENSION(A2D(nn_hls),jpk)    , INTENT(in)   :: pvflx  ! v-flux of advection/diffusion
[7646]320      !
321      INTEGER    ::  ji, jj, jk
[1756]322
[14037]323      IF( cptr /= 'adv' .AND. cptr /= 'ldf' ) RETURN
324      IF( ktra /= jp_tem .AND. ktra /= jp_sal ) RETURN
325
326      IF( cptr == 'adv' ) THEN
327         DO_2D( 0, 0, 0, 0 )
328            hstr_adv(ji,jj,ktra,1) = puflx(ji,jj,1)
329            hstr_adv(ji,jj,ktra,2) = pvflx(ji,jj,1)
330         END_2D
331         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
332            hstr_adv(ji,jj,ktra,1) = hstr_adv(ji,jj,ktra,1) + puflx(ji,jj,jk)
333            hstr_adv(ji,jj,ktra,2) = hstr_adv(ji,jj,ktra,2) + pvflx(ji,jj,jk)
334         END_3D
335      ELSE IF( cptr == 'ldf' ) THEN
336         DO_2D( 0, 0, 0, 0 )
337            hstr_ldf(ji,jj,ktra,1) = puflx(ji,jj,1)
338            hstr_ldf(ji,jj,ktra,2) = pvflx(ji,jj,1)
339         END_2D
340         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
341            hstr_ldf(ji,jj,ktra,1) = hstr_ldf(ji,jj,ktra,1) + puflx(ji,jj,jk)
342            hstr_ldf(ji,jj,ktra,2) = hstr_ldf(ji,jj,ktra,2) + pvflx(ji,jj,jk)
343         END_3D
344      ENDIF
345
346      IF( ntile == 0 .OR. ntile == nijtile ) THEN
347         IF( cptr == 'adv' ) THEN
348            IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,1) )  ! advective heat transport in i-direction
349            IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * hstr_adv(:,:,ktra,1) )  ! advective salt transport in i-direction
350            IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,2) )  ! advective heat transport in j-direction
351            IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * hstr_adv(:,:,ktra,2) )  ! advective salt transport in j-direction
352         ENDIF
353         IF( cptr == 'ldf' ) THEN
354            IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,1) ) ! diffusive heat transport in i-direction
355            IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * hstr_ldf(:,:,ktra,1) ) ! diffusive salt transport in i-direction
356            IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,2) ) ! diffusive heat transport in j-direction
357            IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * hstr_ldf(:,:,ktra,2) ) ! diffusive salt transport in j-direction
358         ENDIF
359      ENDIF
[7646]360   
361   END SUBROUTINE dia_ar5_hst
362
363
[1756]364   SUBROUTINE dia_ar5_init
365      !!----------------------------------------------------------------------
366      !!                  ***  ROUTINE dia_ar5_init  ***
367      !!                   
[2528]368      !! ** Purpose :   initialization for AR5 diagnostic computation
[1756]369      !!----------------------------------------------------------------------
370      INTEGER  ::   inum
[12276]371      INTEGER  ::   ik, idep
[1756]372      INTEGER  ::   ji, jj, jk  ! dummy loop indices
[7753]373      REAL(wp) ::   zztmp 
[9125]374      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
[12276]375      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
[5253]376      !
[1756]377      !!----------------------------------------------------------------------
378      !
[7646]379      l_ar5 = .FALSE.
380      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
381         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
[13089]382         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' ) .OR. &
[14037]383         &  iom_use( 'uadv_heattr' ) .OR. iom_use( 'udiff_heattr' ) .OR. &
384         &  iom_use( 'uadv_salttr' ) .OR. iom_use( 'udiff_salttr' ) .OR. &
385         &  iom_use( 'vadv_heattr' ) .OR. iom_use( 'vdiff_heattr' ) .OR. &
386         &  iom_use( 'vadv_salttr' ) .OR. iom_use( 'vdiff_salttr' ) .OR. &
[13089]387         &  iom_use( 'rhop' )  ) L_ar5 = .TRUE.
[7646]388 
389      IF( l_ar5 ) THEN
390         !
391         !                                      ! allocate dia_ar5 arrays
392         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
[2715]393
[12630]394         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
[1756]395
[12276]396         ALLOCATE( zvol0(jpi,jpj) )
397         zvol0 (:,:) = 0._wp
[7753]398         thick0(:,:) = 0._wp
[14037]399         DO_3D( 1, 1, 1, 1, 1, jpkm1 )   ! interpolation of salinity at the last ocean level (i.e. the partial step)
[12377]400            idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
[12630]401            zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
[12377]402            thick0(ji,jj) = thick0(ji,jj) +  idep   
403         END_3D
[12276]404         vol0 = glob_sum( 'diaar5', zvol0 )
405         DEALLOCATE( zvol0 )
[5253]406
[9125]407         IF( iom_use( 'sshthster' ) ) THEN
[12276]408            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
[9125]409            CALL iom_open ( 'sali_ref_clim_monthly', inum )
[13286]410            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,1), 1  )
411            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,2), 12 )
[9125]412            CALL iom_close( inum )
[6665]413
[9125]414            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
415            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
416            IF( ln_zps ) THEN               ! z-coord. partial steps
[14037]417               DO_2D( 1, 1, 1, 1 )          ! interpolation of salinity at the last ocean level (i.e. the partial step)
[12377]418                  ik = mbkt(ji,jj)
419                  IF( ik > 1 ) THEN
420                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
421                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
422                  ENDIF
423               END_2D
[9125]424            ENDIF
425            !
426            DEALLOCATE( zsaldta )
[7646]427         ENDIF
428         !
[1756]429      ENDIF
430      !
431   END SUBROUTINE dia_ar5_init
432
433   !!======================================================================
434END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.