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 branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 9019

Last change on this file since 9019 was 9019, checked in by timgraham, 6 years ago

Merge of dev_CNRS_2017 into branch

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