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

source: branches/UKMO/dev_r5518_GO6_package_MEDUSA_extra_CMIP6_diags/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 7083

Last change on this file since 7083 was 7083, checked in by malcolmroberts, 8 years ago

Merged in dev_r5518_GO6_package up to revision 7076

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