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

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 7197

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

Manually merge in changes from v3.6_extra_CMIP6_diagnostics branch.
This change also includes a change of the domain_def.xml file so XIOS2 must be used from this revision onwards

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