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

source: branches/UKMO/v3_6_extra_CMIP6_diagnostics/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 7209

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

1) Added in tnpeo changes in diaar5
2) vertical mean temperature fields in field_def.xml

File size: 14.4 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      ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
106      ztsn(:,:,:,jp_sal) = sn0(:,:,:)
107      CALL eos( ztsn, zrhd, fsdept_n(:,:,:) )                       ! now in situ density using initial salinity
108      !
109      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
110      DO jk = 1, jpkm1
111         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
112      END DO
113      IF( .NOT.lk_vvl ) THEN
114         IF ( ln_isfcav ) THEN
115            DO ji=1,jpi
116               DO jj=1,jpj
117                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
118               END DO
119            END DO
120         ELSE
121            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
122         END IF
123      END IF
124      !                                         
125      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
126      IF( lk_mpp )   CALL mpp_sum( zarho )
127      zssh_steric = - zarho / area_tot
128      CALL iom_put( 'sshthster', zssh_steric )
129     
130      !                                         ! steric sea surface height
131      CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) )                 ! now in situ and potential density
132      zrhop(:,:,jpk) = 0._wp
133      CALL iom_put( 'rhop', zrhop )
134      !
135      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
136      DO jk = 1, jpkm1
137         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
138      END DO
139      IF( .NOT.lk_vvl ) THEN
140         IF ( ln_isfcav ) THEN
141            DO ji=1,jpi
142               DO jj=1,jpj
143                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
144               END DO
145            END DO
146         ELSE
147            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
148         END IF
149      END IF
150      !   
151      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
152      IF( lk_mpp )   CALL mpp_sum( zarho )
153      zssh_steric = - zarho / area_tot
154      CALL iom_put( 'sshsteric', zssh_steric )
155     
156      !                                         ! ocean bottom pressure
157      zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
158      zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
159      CALL iom_put( 'botpres', zbotpres )
160
161      !                                         ! Mean density anomalie, temperature and salinity
162      ztemp = 0._wp
163      zsal  = 0._wp
164      DO jk = 1, jpkm1
165         DO jj = 1, jpj
166            DO ji = 1, jpi
167               zztmp = area(ji,jj) * fse3t(ji,jj,jk)
168               ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
169               zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
170            END DO
171         END DO
172      END DO
173      IF( .NOT.lk_vvl ) THEN
174         IF ( ln_isfcav ) THEN
175            DO ji=1,jpi
176               DO jj=1,jpj
177                  ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
178                  zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
179               END DO
180            END DO
181         ELSE
182            ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
183            zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
184         END IF
185      ENDIF
186      IF( lk_mpp ) THEN 
187         CALL mpp_sum( ztemp )
188         CALL mpp_sum( zsal  )
189      END IF
190      !
191      zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
192      ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
193      zsal  = zsal  / zvol                            ! Salinity of liquid seawater
194      !
195      CALL iom_put( 'masstot', zmass )
196      CALL iom_put( 'temptot', ztemp )
197      CALL iom_put( 'saltot' , zsal  )
198
199      IF( iom_use( 'tnpeo' )) THEN   
200      ! Work done against stratification by vertical mixing
201      ! Exclude points where rn2 is negative as convection kicks in here and
202      ! work is not being done against stratification
203          pe(:,:) = 0._wp
204          IF( lk_zdfddm ) THEN
205             DO ji=1,jpi
206                DO jj=1,jpj
207                   DO jk=1,jpk
208                      zrw =   ( fsdepw(ji,jj,jk  ) - fsdept(ji,jj,jk) )   &
209                         &  / ( fsdept(ji,jj,jk-1) - fsdept(ji,jj,jk) )
210                      !
211                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
212                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
213                      !
214                      pe(ji, jj) = pe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * &
215                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
216                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
217
218                   ENDDO
219                ENDDO
220             ENDDO
221          ELSE
222             DO ji=1,jpi
223                DO jj=1,jpj
224                   DO jk=1,jpk
225                       pe(ji,jj) = pe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk)
226                   ENDDO
227                ENDDO
228             ENDDO
229          ENDIF
230          CALL lbc_lnk(pe, 'T', 1._wp)         
231          CALL iom_put( 'tnpeo', pe )
232      ENDIF
233      !
234      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres, pe )
235      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
236      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
237      !
238      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
239      !
240   END SUBROUTINE dia_ar5
241
242
243   SUBROUTINE dia_ar5_init
244      !!----------------------------------------------------------------------
245      !!                  ***  ROUTINE dia_ar5_init  ***
246      !!                   
247      !! ** Purpose :   initialization for AR5 diagnostic computation
248      !!----------------------------------------------------------------------
249      INTEGER  ::   inum
250      INTEGER  ::   ik
251      INTEGER  ::   ji, jj, jk  ! dummy loop indices
252      REAL(wp) ::   zztmp 
253      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
254      ! reading initial file
255      LOGICAL  ::   ln_tsd_init      !: T & S data flag
256      LOGICAL  ::   ln_tsd_tradmp    !: internal damping toward input data flag
257      CHARACTER(len=100)            ::   cn_dir
258      TYPE(FLD_N)                   ::  sn_tem,sn_sal
259      INTEGER  ::   ios=0
260
261      NAMELIST/namtsd/ ln_tsd_init,ln_tsd_tradmp,cn_dir,sn_tem,sn_sal
262      !
263
264      REWIND( numnam_ref )              ! Namelist namtsd in reference namelist :
265      READ  ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901)
266901   IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in reference namelist for dia_ar5', lwp )
267      REWIND( numnam_cfg )              ! Namelist namtsd in configuration namelist : Parameters of the run
268      READ  ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 )
269902   IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in configuration namelist for dia_ar5', lwp )
270      IF(lwm) WRITE ( numond, namtsd )
271      !
272      !!----------------------------------------------------------------------
273      !
274      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
275      !
276      CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
277      !                                      ! allocate dia_ar5 arrays
278      IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
279
280      area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:)
281
282      area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
283
284      vol0        = 0._wp
285      thick0(:,:) = 0._wp
286      DO jk = 1, jpkm1
287         vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
288         thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
289      END DO
290      IF( lk_mpp )   CALL mpp_sum( vol0 )
291
292      CALL iom_open ( TRIM( cn_dir )//TRIM(sn_sal%clname), inum )
293      CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,1), 1  )
294      CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,2), 12 )
295      CALL iom_close( inum )
296      sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
297      sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
298      IF( ln_zps ) THEN               ! z-coord. partial steps
299         DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
300            DO ji = 1, jpi
301               ik = mbkt(ji,jj)
302               IF( ik > 1 ) THEN
303                  zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
304                  sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
305               ENDIF
306            END DO
307         END DO
308      ENDIF
309      !
310      CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
311      !
312      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
313      !
314   END SUBROUTINE dia_ar5_init
315
316#else
317   !!----------------------------------------------------------------------
318   !!   Default option :                                         NO diaar5
319   !!----------------------------------------------------------------------
320   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE.   ! coupled flag
321CONTAINS
322   SUBROUTINE dia_ar5_init    ! Dummy routine
323   END SUBROUTINE dia_ar5_init
324   SUBROUTINE dia_ar5( kt )   ! Empty routine
325      INTEGER ::   kt
326      WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
327   END SUBROUTINE dia_ar5
328#endif
329
330   !!======================================================================
331END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.