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

source: branches/2016/dev_r7233_CMIP6_diags_trunk_version/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 7237

Last change on this file since 7237 was 7237, checked in by timgraham, 8 years ago

Changes in diaar5.F90

  • Property svn:keywords set to Id
File size: 13.3 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#if defined key_diaar5
10   !!----------------------------------------------------------------------
11   !!   'key_diaar5'  :                           activate ar5 diagnotics
12   !!----------------------------------------------------------------------
13   !!   dia_ar5       : AR5 diagnostics
14   !!   dia_ar5_init  : initialisation of AR5 diagnostics
15   !!----------------------------------------------------------------------
16   USE oce            ! ocean dynamics and active tracers
17   USE dom_oce        ! ocean space and time domain
18   USE eosbn2         ! equation of state                (eos_bn2 routine)
19   USE lib_mpp        ! distribued memory computing library
20   USE iom            ! I/O manager library
21   USE timing         ! preformance summary
22   USE wrk_nemo       ! working arrays
23   USE fldread        ! type FLD_N
24   USE phycst         ! physical constant
25   USE in_out_manager  ! I/O manager
26   USE zdfddm
27   USE zdf_oce
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   dia_ar5        ! routine called in step.F90 module
33   PUBLIC   dia_ar5_init   ! routine called in opa.F90 module
34   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
35
36   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE.   ! coupled flag
37
38   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
39   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
40   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   area         ! cell surface (interior domain)
41   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   thick0       ! ocean thickness (interior domain)
42   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity
43     
44   !! * Substitutions
45#  include "zdfddm_substitute.h90"
46   !!----------------------------------------------------------------------
47   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
48   !! $Id$
49   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
50   !!----------------------------------------------------------------------
51CONTAINS
52
53   FUNCTION dia_ar5_alloc()
54      !!----------------------------------------------------------------------
55      !!                    ***  ROUTINE dia_ar5_alloc  ***
56      !!----------------------------------------------------------------------
57      INTEGER :: dia_ar5_alloc
58      !!----------------------------------------------------------------------
59      !
60      ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
61      !
62      IF( lk_mpp             )   CALL mpp_sum ( dia_ar5_alloc )
63      IF( dia_ar5_alloc /= 0 )   CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
64      !
65   END FUNCTION dia_ar5_alloc
66
67
68   SUBROUTINE dia_ar5( kt )
69      !!----------------------------------------------------------------------
70      !!                    ***  ROUTINE dia_ar5  ***
71      !!
72      !! ** Purpose :   compute and output some AR5 diagnostics
73      !!----------------------------------------------------------------------
74      !
75      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
76      !
77      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
78      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
79      !
80      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
81      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace
82      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
83      !!--------------------------------------------------------------------
84      IF( nn_timing == 1 )   CALL timing_start('dia_ar5')
85 
86      CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres, pe )
87      CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    )
88      CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 )
89
90      zarea_ssh(:,:) = area(:,:) * sshn(:,:)
91
92      !                                         ! total volume of liquid seawater
93      zvolssh = SUM( zarea_ssh(:,:) ) 
94      IF( lk_mpp )   CALL mpp_sum( zvolssh )
95      zvol = vol0 + zvolssh
96     
97      CALL iom_put( 'voltot', zvol               )
98      CALL iom_put( 'sshtot', zvolssh / area_tot )
99
100      !                     
101      ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
102      ztsn(:,:,:,jp_sal) = sn0(:,:,:)
103      CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
104      !
105      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
106      DO jk = 1, jpkm1
107         zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
108      END DO
109      IF( ln_linssh ) THEN
110         IF( ln_isfcav ) THEN
111            DO ji=1,jpi
112               DO jj=1,jpj
113                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
114               END DO
115            END DO
116         ELSE
117            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
118         END IF
119!!gm
120!!gm   riceload should be added in both ln_linssh=T or F, no?
121!!gm
122      END IF
123      !                                         
124      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
125      IF( lk_mpp )   CALL mpp_sum( zarho )
126      zssh_steric = - zarho / area_tot
127      CALL iom_put( 'sshthster', zssh_steric )
128     
129      !                                         ! steric sea surface height
130      CALL eos( tsn, zrhd, zrhop, gdept_n(:,:,:) )                 ! now in situ and potential density
131      zrhop(:,:,jpk) = 0._wp
132      CALL iom_put( 'rhop', zrhop )
133      !
134      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
135      DO jk = 1, jpkm1
136         zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk)
137      END DO
138      IF( ln_linssh ) THEN
139         IF ( ln_isfcav ) THEN
140            DO ji=1,jpi
141               DO jj=1,jpj
142                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
143               END DO
144            END DO
145         ELSE
146            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
147         END IF
148      END IF
149      !   
150      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
151      IF( lk_mpp )   CALL mpp_sum( zarho )
152      zssh_steric = - zarho / area_tot
153      CALL iom_put( 'sshsteric', zssh_steric )
154     
155      !                                         ! ocean bottom pressure
156      zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
157      zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
158      CALL iom_put( 'botpres', zbotpres )
159
160      !                                         ! Mean density anomalie, temperature and salinity
161      ztemp = 0._wp
162      zsal  = 0._wp
163      DO jk = 1, jpkm1
164         DO jj = 1, jpj
165            DO ji = 1, jpi
166               zztmp = area(ji,jj) * e3t_n(ji,jj,jk)
167               ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
168               zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
169            END DO
170         END DO
171      END DO
172      IF( ln_linssh ) THEN
173         IF( ln_isfcav ) THEN
174            DO ji=1,jpi
175               DO jj=1,jpj
176                  ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
177                  zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
178               END DO
179            END DO
180         ELSE
181            ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
182            zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
183         END IF
184      ENDIF
185      IF( lk_mpp ) THEN 
186         CALL mpp_sum( ztemp )
187         CALL mpp_sum( zsal  )
188      END IF
189      !
190      zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
191      ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
192      zsal  = zsal  / zvol                            ! Salinity of liquid seawater
193      !
194      CALL iom_put( 'masstot', zmass )
195      CALL iom_put( 'temptot', ztemp )
196      CALL iom_put( 'saltot' , zsal  )
197
198      IF( iom_use( 'tnpeo' )) THEN   
199      ! Work done against stratification by vertical mixing
200      ! Exclude points where rn2 is negative as convection kicks in here and
201      ! work is not being done against stratification
202          pe(:,:) = 0._wp
203          IF( lk_zdfddm ) THEN
204             DO ji=1,jpi
205                DO jj=1,jpj
206                   DO jk=1,jpk
207                      zrw =   ( gdepw_n(ji,jj,jk  ) - gdept_n(ji,jj,jk) )   &
208                         &  / ( gdept_n(ji,jj,jk-1) - gdept_n(ji,jj,jk) )
209                      !
210                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
211                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
212                      !
213                      pe(ji, jj) = pe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * &
214                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
215                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
216
217                   ENDDO
218                ENDDO
219             ENDDO
220          ELSE
221             DO ji=1,jpi
222                DO jj=1,jpj
223                   DO jk=1,jpk
224                       pe(ji,jj) = pe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk)
225                   ENDDO
226                ENDDO
227             ENDDO
228          ENDIF
229          CALL lbc_lnk(pe, 'T', 1._wp)         
230          CALL iom_put( 'tnpeo', pe )
231      ENDIF
232      !
233      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres, pe )
234      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
235      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
236      !
237      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
238      !
239   END SUBROUTINE dia_ar5
240
241
242   SUBROUTINE dia_ar5_init
243      !!----------------------------------------------------------------------
244      !!                  ***  ROUTINE dia_ar5_init  ***
245      !!                   
246      !! ** Purpose :   initialization for AR5 diagnostic computation
247      !!----------------------------------------------------------------------
248      INTEGER  ::   inum
249      INTEGER  ::   ik
250      INTEGER  ::   ji, jj, jk  ! dummy loop indices
251      REAL(wp) ::   zztmp 
252      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
253      !
254      !!----------------------------------------------------------------------
255      !
256      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
257      !
258      CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
259      !                                      ! allocate dia_ar5 arrays
260      IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
261
262      area(:,:) = e1e2t(:,:) * tmask_i(:,:)
263
264      area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
265
266      vol0        = 0._wp
267      thick0(:,:) = 0._wp
268      DO jk = 1, jpkm1
269         vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
270         thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
271      END DO
272      IF( lk_mpp )   CALL mpp_sum( vol0 )
273
274
275      CALL iom_open ( 'sali_ref_clim_monthly', inum )
276      CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
277      CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
278      CALL iom_close( inum )
279
280      sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
281      sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
282      IF( ln_zps ) THEN               ! z-coord. partial steps
283         DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
284            DO ji = 1, jpi
285               ik = mbkt(ji,jj)
286               IF( ik > 1 ) THEN
287                  zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
288                  sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
289               ENDIF
290            END DO
291         END DO
292      ENDIF
293      !
294      CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
295      !
296      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
297      !
298   END SUBROUTINE dia_ar5_init
299
300#else
301   !!----------------------------------------------------------------------
302   !!   Default option :                                         NO diaar5
303   !!----------------------------------------------------------------------
304   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE.   ! coupled flag
305CONTAINS
306   SUBROUTINE dia_ar5_init    ! Dummy routine
307   END SUBROUTINE dia_ar5_init
308   SUBROUTINE dia_ar5( kt )   ! Empty routine
309      INTEGER ::   kt
310      WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
311   END SUBROUTINE dia_ar5
312#endif
313
314   !!======================================================================
315END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.