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 trunk/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 9294

Last change on this file since 9294 was 8083, checked in by timgraham, 7 years ago

Correction to commit r8077 for ticket #1884

  • Property svn:keywords set to Id
File size: 16.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)
[1756]15   USE lib_mpp        ! distribued memory computing library
16   USE iom            ! I/O manager library
[3294]17   USE timing         ! preformance summary
18   USE wrk_nemo       ! working arrays
[5253]19   USE fldread        ! type FLD_N
20   USE phycst         ! physical constant
21   USE in_out_manager  ! I/O manager
[7646]22   USE zdfddm
23   USE zdf_oce
[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 "zdfddm_substitute.h90"
42#  include "vectopt_loop_substitute.h90"
[1756]43   !!----------------------------------------------------------------------
[2528]44   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
45   !! $Id$
46   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[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      !
59      IF( lk_mpp             )   CALL mpp_sum ( dia_ar5_alloc )
60      IF( dia_ar5_alloc /= 0 )   CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
61      !
62   END FUNCTION dia_ar5_alloc
63
64
[1756]65   SUBROUTINE dia_ar5( kt )
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
[2715]73      !
[1756]74      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
75      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
[7646]76      REAL(wp) ::   zaw, zbw, zrw
[3294]77      !
78      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
[7646]79      REAL(wp), POINTER, DIMENSION(:,:)     :: zpe                         ! 2D workspace
[3294]80      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace
81      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
[1756]82      !!--------------------------------------------------------------------
[3294]83      IF( nn_timing == 1 )   CALL timing_start('dia_ar5')
84 
[7646]85      IF( kt == nit000 )     CALL dia_ar5_init
[1756]86
[7646]87      IF( l_ar5 ) THEN
88         CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres )
89         CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    )
90         CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 )
[7753]91         zarea_ssh(:,:) = area(:,:) * sshn(:,:)
[7646]92      ENDIF
93      !
94      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
95         !                                         ! total volume of liquid seawater
96         zvolssh = SUM( zarea_ssh(:,:) ) 
97         IF( lk_mpp )   CALL mpp_sum( zvolssh )
98         zvol = vol0 + zvolssh
[1756]99     
[7646]100         CALL iom_put( 'voltot', zvol               )
101         CALL iom_put( 'sshtot', zvolssh / area_tot )
102         CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
103         !
104      ENDIF
[1756]105
[7646]106      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN   
107         !                     
[7753]108         ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
109         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
[7646]110         CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
111         !
[7753]112         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]113         DO jk = 1, jpkm1
[7753]114            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
[7646]115         END DO
116         IF( ln_linssh ) THEN
117            IF( ln_isfcav ) THEN
118               DO ji = 1, jpi
119                  DO jj = 1, jpj
120                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
121                  END DO
[5120]122               END DO
[7646]123            ELSE
[7753]124               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
[7646]125            END IF
[6140]126!!gm
127!!gm   riceload should be added in both ln_linssh=T or F, no?
128!!gm
[7646]129         END IF
130         !                                         
[7753]131         zarho = SUM( area(:,:) * zbotpres(:,:) ) 
[7646]132         IF( lk_mpp )   CALL mpp_sum( zarho )
133         zssh_steric = - zarho / area_tot
134         CALL iom_put( 'sshthster', zssh_steric )
[1756]135     
[7646]136         !                                         ! steric sea surface height
137         CALL eos( tsn, zrhd, zrhop, gdept_n(:,:,:) )                 ! now in situ and potential density
[7753]138         zrhop(:,:,jpk) = 0._wp
[7646]139         CALL iom_put( 'rhop', zrhop )
140         !
[7753]141         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
[7646]142         DO jk = 1, jpkm1
[7753]143            zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
[7646]144         END DO
145         IF( ln_linssh ) THEN
146            IF ( ln_isfcav ) THEN
147               DO ji = 1,jpi
148                  DO jj = 1,jpj
149                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
150                  END DO
[5120]151               END DO
[7646]152            ELSE
[7753]153               zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
[7646]154            END IF
[5120]155         END IF
[7646]156         !   
[7753]157         zarho = SUM( area(:,:) * zbotpres(:,:) ) 
[7646]158         IF( lk_mpp )   CALL mpp_sum( zarho )
159         zssh_steric = - zarho / area_tot
160         CALL iom_put( 'sshsteric', zssh_steric )
[1756]161     
[7646]162         !                                         ! ocean bottom pressure
163         zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
[7753]164         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
[7646]165         CALL iom_put( 'botpres', zbotpres )
166         !
167      ENDIF
[1756]168
[7646]169      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
170         !                                         ! Mean density anomalie, temperature and salinity
171         ztemp = 0._wp
172         zsal  = 0._wp
173         DO jk = 1, jpkm1
174            DO jj = 1, jpj
175               DO ji = 1, jpi
176                  zztmp = area(ji,jj) * e3t_n(ji,jj,jk)
177                  ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
178                  zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
179               END DO
[1756]180            END DO
181         END DO
[7646]182         IF( ln_linssh ) THEN
183            IF( ln_isfcav ) THEN
184               DO ji = 1, jpi
185                  DO jj = 1, jpj
186                     ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
187                     zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
188                  END DO
[5120]189               END DO
[7646]190            ELSE
191               ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
192               zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
193            END IF
194         ENDIF
195         IF( lk_mpp ) THEN 
196            CALL mpp_sum( ztemp )
197            CALL mpp_sum( zsal  )
[5120]198         END IF
[7646]199         !
200         zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
201         ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
202         zsal  = zsal  / zvol                            ! Salinity of liquid seawater
203         !
204         CALL iom_put( 'masstot', zmass )
205         CALL iom_put( 'temptot', ztemp )
206         CALL iom_put( 'saltot' , zsal  )
207         !
[1756]208      ENDIF
[7646]209
210      IF( iom_use( 'tnpeo' )) THEN   
211      ! Work done against stratification by vertical mixing
212      ! Exclude points where rn2 is negative as convection kicks in here and
213      ! work is not being done against stratification
[8078]214         CALL wrk_alloc( jpi, jpj, zpe )
215         zpe(:,:) = 0._wp
[8083]216         IF( lk_zdfddm ) THEN
[8078]217            DO jk = 2, jpk
218               DO jj = 1, jpj
219                  DO ji = 1, jpi
220                     IF( rn2(ji,jj,jk) > 0._wp ) THEN
221                        zrw =   ( gdepw_n(ji,jj,jk  ) - gdept_n(ji,jj,jk) )   &
222                           &  / ( gdept_n(ji,jj,jk-1) - gdept_n(ji,jj,jk) )
223                        !
224                        zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
225                        zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
226                        !
227                        zpe(ji, jj) = zpe(ji, jj)            &
228                           &        -  grav * (    avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
229                           &                   - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
230                     ENDIF
231                  END DO
232               END DO
233             END DO
[7646]234          ELSE
[8078]235            DO jk = 1, jpk
236               DO ji = 1, jpi
237                  DO jj = 1, jpj
238                     zpe(ji,jj) = zpe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * e3w_n(ji, jj, jk)
239                  END DO
240               END DO
241            END DO
242         ENDIF
243         CALL lbc_lnk( zpe, 'T', 1._wp)         
244         CALL iom_put( 'tnpeo', zpe )
245         CALL wrk_dealloc( jpi, jpj, zpe )
[7646]246      ENDIF
[1756]247      !
[7646]248      IF( l_ar5 ) THEN
249        CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres )
250        CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
251        CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
252      ENDIF
[1756]253      !
[3294]254      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
255      !
[1756]256   END SUBROUTINE dia_ar5
257
[7646]258   SUBROUTINE dia_ar5_hst( ktra, cptr, pua, pva ) 
259      !!----------------------------------------------------------------------
260      !!                    ***  ROUTINE dia_ar5_htr ***
261      !!----------------------------------------------------------------------
262      !! Wrapper for heat transport calculations
263      !! Called from all advection and/or diffusion routines
264      !!----------------------------------------------------------------------
265      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
266      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
267      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pua   ! 3D input array of advection/diffusion
268      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pva   ! 3D input array of advection/diffusion
269      !
270      INTEGER    ::  ji, jj, jk
271      REAL(wp), POINTER, DIMENSION(:,:)  :: z2d
[1756]272
[7646]273   
274
275      CALL wrk_alloc( jpi, jpj, z2d )
276      z2d(:,:) = pua(:,:,1) 
277      DO jk = 1, jpkm1
278         DO jj = 2, jpjm1
279            DO ji = fs_2, fs_jpim1   ! vector opt.
280               z2d(ji,jj) = z2d(ji,jj) + pua(ji,jj,jk) 
281            END DO
282         END DO
283       END DO
284       CALL lbc_lnk( z2d, 'U', -1. )
285       IF( cptr == 'adv' ) THEN
286          IF( ktra == jp_tem ) CALL iom_put( "uadv_heattr" , rau0_rcp * z2d )  ! advective heat transport in i-direction
287          IF( ktra == jp_sal ) CALL iom_put( "uadv_salttr" , rau0     * z2d )  ! advective salt transport in i-direction
288       ENDIF
289       IF( cptr == 'ldf' ) THEN
290          IF( ktra == jp_tem ) CALL iom_put( "udiff_heattr" , rau0_rcp * z2d ) ! diffusive heat transport in i-direction
291          IF( ktra == jp_sal ) CALL iom_put( "udiff_salttr" , rau0     * z2d ) ! diffusive salt transport in i-direction
292       ENDIF
293       !
294       z2d(:,:) = pva(:,:,1) 
295       DO jk = 1, jpkm1
296          DO jj = 2, jpjm1
297             DO ji = fs_2, fs_jpim1   ! vector opt.
298                z2d(ji,jj) = z2d(ji,jj) + pva(ji,jj,jk) 
299             END DO
300          END DO
301       END DO
302       CALL lbc_lnk( z2d, 'V', -1. )
303       IF( cptr == 'adv' ) THEN
304          IF( ktra == jp_tem ) CALL iom_put( "vadv_heattr" , rau0_rcp * z2d )  ! advective heat transport in j-direction
305          IF( ktra == jp_sal ) CALL iom_put( "vadv_salttr" , rau0     * z2d )  ! advective salt transport in j-direction
306       ENDIF
307       IF( cptr == 'ldf' ) THEN
308          IF( ktra == jp_tem ) CALL iom_put( "vdiff_heattr" , rau0_rcp * z2d ) ! diffusive heat transport in j-direction
309          IF( ktra == jp_sal ) CALL iom_put( "vdiff_salttr" , rau0     * z2d ) ! diffusive salt transport in j-direction
310       ENDIF
311         
312       CALL wrk_dealloc( jpi, jpj, z2d )
313
314   END SUBROUTINE dia_ar5_hst
315
316
[1756]317   SUBROUTINE dia_ar5_init
318      !!----------------------------------------------------------------------
319      !!                  ***  ROUTINE dia_ar5_init  ***
320      !!                   
[2528]321      !! ** Purpose :   initialization for AR5 diagnostic computation
[1756]322      !!----------------------------------------------------------------------
323      INTEGER  ::   inum
324      INTEGER  ::   ik
325      INTEGER  ::   ji, jj, jk  ! dummy loop indices
[7753]326      REAL(wp) ::   zztmp 
[2715]327      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
[5253]328      !
[1756]329      !!----------------------------------------------------------------------
330      !
[3294]331      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
332      !
[7646]333      l_ar5 = .FALSE.
334      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
335         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
336         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) L_ar5 = .TRUE.
337 
338      IF( l_ar5 ) THEN
339         !
340         CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
341         !                                      ! allocate dia_ar5 arrays
342         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
[2715]343
[7753]344         area(:,:) = e1e2t(:,:) * tmask_i(:,:)
[1756]345
[7646]346         area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
[1756]347
[7646]348         vol0        = 0._wp
[7753]349         thick0(:,:) = 0._wp
[7646]350         DO jk = 1, jpkm1
[7753]351            vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
352            thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
[7646]353         END DO
354         IF( lk_mpp )   CALL mpp_sum( vol0 )
[5253]355
[7646]356         CALL iom_open ( 'sali_ref_clim_monthly', inum )
357         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
358         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
359         CALL iom_close( inum )
[6665]360
[7753]361         sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
362         sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
[7646]363         IF( ln_zps ) THEN               ! z-coord. partial steps
364            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
365               DO ji = 1, jpi
366                  ik = mbkt(ji,jj)
367                  IF( ik > 1 ) THEN
368                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
369                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
370                  ENDIF
371               END DO
[1756]372            END DO
[7646]373         ENDIF
374         !
375         CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
376         !
[1756]377      ENDIF
378      !
[3294]379      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
380      !
[1756]381   END SUBROUTINE dia_ar5_init
382
383   !!======================================================================
384END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.