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

source: branches/UKMO/dev_r5518_GO6_under_ice_relax_dr_hook/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 11738

Last change on this file since 11738 was 11738, checked in by marc, 5 years ago

The Dr Hook changes from my perl code.

File size: 15.6 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   USE yomhook, ONLY: lhook, dr_hook
30   USE parkind1, ONLY: jprb, jpim
31
32   IMPLICIT NONE
33   PRIVATE
34
35   PUBLIC   dia_ar5        ! routine called in step.F90 module
36   PUBLIC   dia_ar5_init   ! routine called in opa.F90 module
37   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
38
39   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE.   ! coupled flag
40
41   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
42   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
43   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   area         ! cell surface (interior domain)
44   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   thick0       ! ocean thickness (interior domain)
45   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity
46     
47   !! * Substitutions
48#  include "domzgr_substitute.h90"
49#  include "zdfddm_substitute.h90"
50   !!----------------------------------------------------------------------
51   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
52   !! $Id$
53   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
54   !!----------------------------------------------------------------------
55CONTAINS
56
57   FUNCTION dia_ar5_alloc()
58      !!----------------------------------------------------------------------
59      !!                    ***  ROUTINE dia_ar5_alloc  ***
60      !!----------------------------------------------------------------------
61      INTEGER :: dia_ar5_alloc
62      INTEGER(KIND=jpim), PARAMETER :: zhook_in = 0
63      INTEGER(KIND=jpim), PARAMETER :: zhook_out = 1
64      REAL(KIND=jprb)               :: zhook_handle
65
66      CHARACTER(LEN=*), PARAMETER :: RoutineName='DIA_AR5_ALLOC'
67
68      IF (lhook) CALL dr_hook(RoutineName,zhook_in,zhook_handle)
69
70      !!----------------------------------------------------------------------
71      !
72      ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
73      !
74      IF( lk_mpp             )   CALL mpp_sum ( dia_ar5_alloc )
75      IF( dia_ar5_alloc /= 0 )   CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
76      !
77      IF (lhook) CALL dr_hook(RoutineName,zhook_out,zhook_handle)
78   END FUNCTION dia_ar5_alloc
79
80
81   SUBROUTINE dia_ar5( kt )
82      !!----------------------------------------------------------------------
83      !!                    ***  ROUTINE dia_ar5  ***
84      !!
85      !! ** Purpose :   compute and output some AR5 diagnostics
86      !!----------------------------------------------------------------------
87      !
88      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
89      !
90      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
91      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
92      REAL(wp) ::   zaw, zbw, zrw
93      !
94      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
95      REAL(wp), POINTER, DIMENSION(:,:)     :: pe                         ! 2D workspace
96      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace
97      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
98      INTEGER(KIND=jpim), PARAMETER :: zhook_in = 0
99      INTEGER(KIND=jpim), PARAMETER :: zhook_out = 1
100      REAL(KIND=jprb)               :: zhook_handle
101
102      CHARACTER(LEN=*), PARAMETER :: RoutineName='DIA_AR5'
103
104      IF (lhook) CALL dr_hook(RoutineName,zhook_in,zhook_handle)
105
106      !!--------------------------------------------------------------------
107      IF( nn_timing == 1 )   CALL timing_start('dia_ar5')
108
109      !Call to init moved to here so that we can call iom_use in the
110      !initialisation
111      IF( kt == nit000 )     CALL dia_ar5_init
112 
113      CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres, pe )
114      CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    )
115      CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 )
116
117      zarea_ssh(:,:) = area(:,:) * sshn(:,:)
118
119      !                                         ! total volume of liquid seawater
120      zvolssh = SUM( zarea_ssh(:,:) ) 
121      IF( lk_mpp )   CALL mpp_sum( zvolssh )
122      zvol = vol0 + zvolssh
123     
124      CALL iom_put( 'voltot', zvol               )
125      CALL iom_put( 'sshtot', zvolssh / area_tot )
126      CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
127
128      !                     
129      IF( iom_use('sshthster')) THEN
130         ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
131         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
132         CALL eos( ztsn, zrhd, fsdept_n(:,:,:) )                       ! now in situ density using initial salinity
133         !
134         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
135         DO jk = 1, jpkm1
136            zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
137         END DO
138         IF( .NOT.lk_vvl ) 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( 'sshthster', zssh_steric )
154      ENDIF
155     
156      !                                         ! steric sea surface height
157      CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) )                 ! now in situ and potential density
158      zrhop(:,:,jpk) = 0._wp
159      CALL iom_put( 'rhop', zrhop )
160      !
161      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
162      DO jk = 1, jpkm1
163         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
164      END DO
165      IF( .NOT.lk_vvl ) THEN
166         IF ( ln_isfcav ) THEN
167            DO ji=1,jpi
168               DO jj=1,jpj
169                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
170               END DO
171            END DO
172         ELSE
173            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
174         END IF
175      END IF
176      !   
177      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
178      IF( lk_mpp )   CALL mpp_sum( zarho )
179      zssh_steric = - zarho / area_tot
180      CALL iom_put( 'sshsteric', zssh_steric )
181     
182      !                                         ! ocean bottom pressure
183      zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
184      zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
185      CALL iom_put( 'botpres', zbotpres )
186
187      !                                         ! Mean density anomalie, temperature and salinity
188      ztemp = 0._wp
189      zsal  = 0._wp
190      DO jk = 1, jpkm1
191         DO jj = 1, jpj
192            DO ji = 1, jpi
193               zztmp = area(ji,jj) * fse3t(ji,jj,jk)
194               ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
195               zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
196            END DO
197         END DO
198      END DO
199      IF( .NOT.lk_vvl ) THEN
200         IF ( ln_isfcav ) THEN
201            DO ji=1,jpi
202               DO jj=1,jpj
203                  ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
204                  zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
205               END DO
206            END DO
207         ELSE
208            ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
209            zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
210         END IF
211      ENDIF
212      IF( lk_mpp ) THEN 
213         CALL mpp_sum( ztemp )
214         CALL mpp_sum( zsal  )
215      END IF
216      !
217      zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
218      ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
219      zsal  = zsal  / zvol                            ! Salinity of liquid seawater
220      !
221      CALL iom_put( 'masstot', zmass )
222      CALL iom_put( 'temptot', ztemp )
223      CALL iom_put( 'saltot' , zsal  )
224
225      IF( iom_use( 'tnpeo' )) THEN   
226      ! Work done against stratification by vertical mixing
227      ! Exclude points where rn2 is negative as convection kicks in here and
228      ! work is not being done against stratification
229          pe(:,:) = 0._wp
230          IF( lk_zdfddm ) THEN
231             DO ji=1,jpi
232                DO jj=1,jpj
233                   DO jk=1,jpk
234                      zrw =   ( fsdepw(ji,jj,jk  ) - fsdept(ji,jj,jk) )   &
235                         &  / ( fsdept(ji,jj,jk-1) - fsdept(ji,jj,jk) )
236                      !
237                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
238                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
239                      !
240                      pe(ji, jj) = pe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * &
241                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
242                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
243
244                   ENDDO
245                ENDDO
246             ENDDO
247          ELSE
248             DO ji=1,jpi
249                DO jj=1,jpj
250                   DO jk=1,jpk
251                       pe(ji,jj) = pe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk)
252                   ENDDO
253                ENDDO
254             ENDDO
255          ENDIF
256          CALL iom_put( 'tnpeo', pe )
257      ENDIF
258      !
259      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres, pe )
260      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
261      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
262      !
263      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
264      !
265      IF (lhook) CALL dr_hook(RoutineName,zhook_out,zhook_handle)
266   END SUBROUTINE dia_ar5
267
268
269   SUBROUTINE dia_ar5_init
270      !!----------------------------------------------------------------------
271      !!                  ***  ROUTINE dia_ar5_init  ***
272      !!                   
273      !! ** Purpose :   initialization for AR5 diagnostic computation
274      !!----------------------------------------------------------------------
275      INTEGER  ::   inum
276      INTEGER  ::   ik
277      INTEGER  ::   ji, jj, jk  ! dummy loop indices
278      REAL(wp) ::   zztmp 
279      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
280      INTEGER(KIND=jpim), PARAMETER :: zhook_in = 0
281      INTEGER(KIND=jpim), PARAMETER :: zhook_out = 1
282      REAL(KIND=jprb)               :: zhook_handle
283
284      CHARACTER(LEN=*), PARAMETER :: RoutineName='DIA_AR5_INIT'
285
286      IF (lhook) CALL dr_hook(RoutineName,zhook_in,zhook_handle)
287
288      !
289      !!----------------------------------------------------------------------
290      !
291      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
292      !
293      CALL wrk_alloc( jpi, jpj, jpk, 2, zsaldta )
294      !                                      ! allocate dia_ar5 arrays
295      IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
296
297      area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:)
298
299      area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
300
301      vol0        = 0._wp
302      thick0(:,:) = 0._wp
303      DO jk = 1, jpkm1
304         vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
305         thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
306      END DO
307      IF( lk_mpp )   CALL mpp_sum( vol0 )
308
309      IF( iom_use('sshthster')) THEN
310         CALL iom_open ( 'sali_ref_clim_monthly', inum )
311         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
312         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
313         CALL iom_close( inum )
314
315         sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
316         sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
317         IF( ln_zps ) THEN               ! z-coord. partial steps
318            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
319               DO ji = 1, jpi
320                  ik = mbkt(ji,jj)
321                  IF( ik > 1 ) THEN
322                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
323                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
324                  ENDIF
325               END DO
326            END DO
327         ENDIF
328      ENDIF
329      !
330      CALL wrk_dealloc( jpi, jpj, jpk, 2, zsaldta )
331      !
332      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
333      !
334      IF (lhook) CALL dr_hook(RoutineName,zhook_out,zhook_handle)
335   END SUBROUTINE dia_ar5_init
336
337#else
338   !!----------------------------------------------------------------------
339   !!   Default option :                                         NO diaar5
340   !!----------------------------------------------------------------------
341   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE.   ! coupled flag
342CONTAINS
343   SUBROUTINE dia_ar5_init    ! Dummy routine
344   INTEGER(KIND=jpim), PARAMETER :: zhook_in = 0
345   INTEGER(KIND=jpim), PARAMETER :: zhook_out = 1
346   REAL(KIND=jprb)               :: zhook_handle
347
348   CHARACTER(LEN=*), PARAMETER :: RoutineName='DIA_AR5_INIT'
349
350   IF (lhook) CALL dr_hook(RoutineName,zhook_in,zhook_handle)
351
352   IF (lhook) CALL dr_hook(RoutineName,zhook_out,zhook_handle)
353   END SUBROUTINE dia_ar5_init
354   SUBROUTINE dia_ar5( kt )   ! Empty routine
355      INTEGER ::   kt
356      INTEGER(KIND=jpim), PARAMETER :: zhook_in = 0
357      INTEGER(KIND=jpim), PARAMETER :: zhook_out = 1
358      REAL(KIND=jprb)               :: zhook_handle
359
360      CHARACTER(LEN=*), PARAMETER :: RoutineName='DIA_AR5'
361
362      IF (lhook) CALL dr_hook(RoutineName,zhook_in,zhook_handle)
363
364      WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
365      IF (lhook) CALL dr_hook(RoutineName,zhook_out,zhook_handle)
366   END SUBROUTINE dia_ar5
367#endif
368
369   !!======================================================================
370END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.