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

source: branches/UKMO/AMM15_v3_6_STABLE_package_reanalysis2/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 10962

Last change on this file since 10962 was 10962, checked in by rrenshaw, 5 years ago

code fix

File size: 15.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
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   dia_ar5        ! routine called in step.F90 module
31   PUBLIC   dia_ar5_init   ! routine called in opa.F90 module
32   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
33
34   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE.   ! coupled flag
35
36   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
37   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
38   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   area         ! cell surface (interior domain)
39   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   thick0       ! ocean thickness (interior domain)
40   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity
41   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   tn0          ! initial temperature
42   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshthster_mat         ! ssh_thermosteric height
43   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshhlster_mat         ! ssh_halosteric height
44   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshsteric_mat         ! ssh_steric height
45   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   zbotpres_mat          ! bottom pressure
46     
47   !! * Substitutions
48#  include "domzgr_substitute.h90"
49   !!----------------------------------------------------------------------
50   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
51   !! $Id$
52   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
53   !!----------------------------------------------------------------------
54CONTAINS
55
56   FUNCTION dia_ar5_alloc()
57      !!----------------------------------------------------------------------
58      !!                    ***  ROUTINE dia_ar5_alloc  ***
59      !!----------------------------------------------------------------------
60      INTEGER :: dia_ar5_alloc
61      !!----------------------------------------------------------------------
62      !
63      ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk), tn0(jpi,jpj,jpk) , &
64          & sshthster_mat(jpi,jpj),sshhlster_mat(jpi,jpj),sshsteric_mat(jpi,jpj), &
65          & zbotpres_mat(jpi,jpj),STAT=dia_ar5_alloc )
66      !
67      IF( lk_mpp             )   CALL mpp_sum ( dia_ar5_alloc )
68      IF( dia_ar5_alloc /= 0 )   CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
69      !
70   END FUNCTION dia_ar5_alloc
71
72
73   SUBROUTINE dia_ar5( kt )
74      !!----------------------------------------------------------------------
75      !!                    ***  ROUTINE dia_ar5  ***
76      !!
77      !! ** Purpose :   compute and output some AR5 diagnostics
78      !!----------------------------------------------------------------------
79      !
80      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
81      !
82      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
83      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
84      !
85      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
86      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace
87      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
88      !!--------------------------------------------------------------------
89      IF( nn_timing == 1 )   CALL timing_start('dia_ar5')
90 
91      CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres )
92      CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    )
93      CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 )
94     
95      sshthster_mat(:,:) = 0._wp 
96      sshhlster_mat(:,:) = 0._wp 
97      sshsteric_mat(:,:) = 0._wp 
98      zbotpres_mat(:,:)  = 0._wp 
99
100      zarea_ssh(:,:) = area(:,:) * sshn(:,:)
101
102      !                                         ! total volume of liquid seawater
103      zvolssh = SUM( zarea_ssh(:,:) ) 
104      IF( lk_mpp )   CALL mpp_sum( zvolssh )
105      zvol = vol0 + zvolssh
106     
107      CALL iom_put( 'voltot', zvol               )
108      CALL iom_put( 'sshtot', zvolssh / area_tot )
109
110      !                     
111      ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
112      ztsn(:,:,:,jp_sal) = sn0(:,:,:)
113      CALL eos( ztsn, zrhd, fsdept_n(:,:,:) )                       ! now in situ density using initial salinity
114      !
115      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
116      DO jk = 1, jpkm1
117         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
118      END DO
119      IF( .NOT.lk_vvl ) THEN
120         IF ( ln_isfcav ) THEN
121            DO ji=1,jpi
122               DO jj=1,jpj
123                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
124               END DO
125            END DO
126         ELSE
127            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
128         END IF
129      END IF
130      !                                         
131      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
132      IF( lk_mpp )   CALL mpp_sum( zarho )
133      zssh_steric = - zarho / area_tot
134      CALL iom_put( 'sshthster', zssh_steric )
135      sshthster_mat(:,:) =  -zbotpres(:,:)
136      CALL iom_put( 'sshthster_mat', sshthster_mat )
137
138      !                     
139      ztsn(:,:,:,jp_tem) = tn0(:,:,:)                    ! thermohaline ssh
140      ztsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) 
141      CALL eos( ztsn, zrhd, fsdept_n(:,:,:) )                       ! now in situ density using initial temperature
142      !
143      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
144      DO jk = 1, jpkm1
145         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
146      END DO
147      IF( .NOT.lk_vvl ) THEN
148         IF ( ln_isfcav ) THEN
149            DO ji=1,jpi
150               DO jj=1,jpj
151                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
152               END DO
153            END DO
154         ELSE
155            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
156         END IF
157      END IF
158      !                                         
159      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
160      IF( lk_mpp )   CALL mpp_sum( zarho )
161      zssh_steric = - zarho / area_tot
162      CALL iom_put( 'sshhlster', zssh_steric )
163      sshhlster_mat(:,:) = -zbotpres(:,:)
164      CALL iom_put( 'sshhlster_mat', sshhlster_mat )
165     
166
167     
168      !                                         ! steric sea surface height
169      CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) )                 ! now in situ and potential density
170      zrhop(:,:,jpk) = 0._wp
171      CALL iom_put( 'rhop', zrhop )
172      !
173      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
174      DO jk = 1, jpkm1
175         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
176      END DO
177      IF( .NOT.lk_vvl ) THEN
178         IF ( ln_isfcav ) THEN
179            DO ji=1,jpi
180               DO jj=1,jpj
181                  zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
182               END DO
183            END DO
184         ELSE
185            zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
186         END IF
187      END IF
188      !   
189      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
190      IF( lk_mpp )   CALL mpp_sum( zarho )
191      zssh_steric = - zarho / area_tot
192      CALL iom_put( 'sshsteric', zssh_steric )
193      sshsteric_mat(:,:) = -zbotpres(:,:) 
194      CALL iom_put( 'sshsteric_mat', sshsteric_mat )
195     
196      !                                         ! ocean bottom pressure
197      zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
198      zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
199      zbotpres_mat(:,:) = zbotpres(:,:)
200      CALL iom_put( 'botpres', zbotpres_mat )
201     
202      !                                         ! Mean density anomalie, temperature and salinity
203      ztemp = 0._wp
204      zsal  = 0._wp
205      DO jk = 1, jpkm1
206         DO jj = 1, jpj
207            DO ji = 1, jpi
208               zztmp = area(ji,jj) * fse3t(ji,jj,jk)
209               ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
210               zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
211            END DO
212         END DO
213      END DO
214      IF( .NOT.lk_vvl ) THEN
215         IF ( ln_isfcav ) THEN
216            DO ji=1,jpi
217               DO jj=1,jpj
218                  ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
219                  zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
220               END DO
221            END DO
222         ELSE
223            ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
224            zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
225         END IF
226      ENDIF
227      IF( lk_mpp ) THEN 
228         CALL mpp_sum( ztemp )
229         CALL mpp_sum( zsal  )
230      END IF
231      !
232      zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
233      ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
234      zsal  = zsal  / zvol                            ! Salinity of liquid seawater
235      !
236      CALL iom_put( 'masstot', zmass )
237      CALL iom_put( 'temptot', ztemp )
238      CALL iom_put( 'saltot' , zsal  )
239      !
240      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres )
241      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
242      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
243      !
244      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
245      !
246   END SUBROUTINE dia_ar5
247
248
249   SUBROUTINE dia_ar5_init
250      !!----------------------------------------------------------------------
251      !!                  ***  ROUTINE dia_ar5_init  ***
252      !!                   
253      !! ** Purpose :   initialization for AR5 diagnostic computation
254      !!----------------------------------------------------------------------
255      INTEGER  ::   inum
256      INTEGER  ::   ik
257      INTEGER  ::   ji, jj, jk  ! dummy loop indices
258      REAL(wp) ::   zztmp 
259      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
260      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   ztemdta   ! Jan/Dec levitus salinity
261      ! reading initial file
262      LOGICAL  ::   ln_tsd_init      !: T & S data flag
263      LOGICAL  ::   ln_tsd_tradmp    !: internal damping toward input data flag
264      CHARACTER(len=100)            ::   cn_dir
265      TYPE(FLD_N)                   ::  sn_tem,sn_sal
266      INTEGER  ::   ios=0
267
268      NAMELIST/namtsd/ ln_tsd_init,ln_tsd_tradmp,cn_dir,sn_tem,sn_sal
269      !
270
271      REWIND( numnam_ref )              ! Namelist namtsd in reference namelist :
272      READ  ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901)
273901   IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in reference namelist for dia_ar5', lwp )
274      REWIND( numnam_cfg )              ! Namelist namtsd in configuration namelist : Parameters of the run
275      READ  ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 )
276902   IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in configuration namelist for dia_ar5', lwp )
277      IF(lwm) WRITE ( numond, namtsd )
278      !
279      !!----------------------------------------------------------------------
280      !
281      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
282      !
283      CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
284      CALL wrk_alloc( jpi , jpj , jpk, jpts, ztemdta )
285      !                                      ! allocate dia_ar5 arrays
286      IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
287
288      area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:)
289
290      area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
291
292      vol0        = 0._wp
293      thick0(:,:) = 0._wp
294      DO jk = 1, jpkm1
295         vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
296         thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
297      END DO
298      IF( lk_mpp )   CALL mpp_sum( vol0 )
299
300      CALL iom_open ( TRIM( cn_dir )//TRIM(sn_sal%clname), inum )
301      CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,1), 1  )
302      CALL iom_get  ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,2), 12 )
303      CALL iom_close( inum )
304
305      CALL iom_open ( TRIM( cn_dir )//TRIM(sn_tem%clname), inum )
306      CALL iom_get  ( inum, jpdom_data, TRIM(sn_tem%clvar), ztemdta(:,:,:,1), 1  )
307      CALL iom_get  ( inum, jpdom_data, TRIM(sn_tem%clvar), ztemdta(:,:,:,2), 12 )
308      CALL iom_close( inum )
309      sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
310      tn0(:,:,:) = 0.5_wp * ( ztemdta(:,:,:,1) + ztemdta(:,:,:,2) )       
311      sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) 
312      tn0(:,:,:) = tn0(:,:,:) * tmask(:,:,:)
313      IF( ln_zps ) THEN               ! z-coord. partial steps
314         DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
315            DO ji = 1, jpi
316               ik = mbkt(ji,jj)
317               IF( ik > 1 ) THEN
318                  zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
319                  sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
320                  tn0(ji,jj,ik) = ( 1._wp - zztmp ) * tn0(ji,jj,ik) + zztmp * tn0(ji,jj,ik-1)
321               ENDIF
322            END DO
323         END DO
324      ENDIF
325      !
326      CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
327      CALL wrk_dealloc( jpi , jpj , jpk, jpts, ztemdta )
328      !
329      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
330      !
331   END SUBROUTINE dia_ar5_init
332
333#else
334   !!----------------------------------------------------------------------
335   !!   Default option :                                         NO diaar5
336   !!----------------------------------------------------------------------
337   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE.   ! coupled flag
338CONTAINS
339   SUBROUTINE dia_ar5_init    ! Dummy routine
340   END SUBROUTINE dia_ar5_init
341   SUBROUTINE dia_ar5( kt )   ! Empty routine
342      INTEGER ::   kt
343      WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
344   END SUBROUTINE dia_ar5
345#endif
346
347   !!======================================================================
348END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.