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

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 5866

Last change on this file since 5866 was 5866, checked in by gm, 8 years ago

#1613: vvl by default: add ln_linssh and remove key_vvl

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