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

source: trunk/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 5107

Last change on this file since 5107 was 5107, checked in by smasson, 9 years ago

update IO for VVL compatibility, see ticket #1474

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