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

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 4292

Last change on this file since 4292 was 4292, checked in by cetlod, 10 years ago

dev_MERGE_2013 : 1st step of the merge, see ticket #1185

  • Property svn:keywords set to Id
File size: 10.7 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      CALL iom_put( 'cellthc', fse3t(:,:,:) )
86
87      zarea_ssh(:,:) = area(:,:) * sshn(:,:)
88
89      !                                         ! total volume of liquid seawater
90      zvolssh = SUM( zarea_ssh(:,:) ) 
91      IF( lk_mpp )   CALL mpp_sum( zvolssh )
92      zvol = vol0 + zvolssh
93     
94      CALL iom_put( 'voltot', zvol               )
95      CALL iom_put( 'sshtot', zvolssh / area_tot )
96
97      !                     
98      ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)                    ! thermosteric ssh
99      ztsn(:,:,:,jp_sal) = sn0(:,:,:)
100      CALL eos( ztsn, zrhd )                       ! now in situ density using initial salinity
101      !
102      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
103      DO jk = 1, jpkm1
104         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
105      END DO
106      IF( .NOT.lk_vvl )   zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
107      !                                         
108      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
109      IF( lk_mpp )   CALL mpp_sum( zarho )
110      zssh_steric = - zarho / area_tot
111      CALL iom_put( 'sshthster', zssh_steric )
112     
113      !                                         ! steric sea surface height
114      CALL eos( tsn, zrhd, zrhop )                 ! now in situ and potential density
115      zrhop(:,:,jpk) = 0._wp
116      CALL iom_put( 'rhop', zrhop )
117      !
118      zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
119      DO jk = 1, jpkm1
120         zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
121      END DO
122      IF( .NOT.lk_vvl )   zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
123      !   
124      zarho = SUM( area(:,:) * zbotpres(:,:) ) 
125      IF( lk_mpp )   CALL mpp_sum( zarho )
126      zssh_steric = - zarho / area_tot
127      CALL iom_put( 'sshsteric', zssh_steric )
128     
129      !                                         ! ocean bottom pressure
130      zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
131      zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
132      CALL iom_put( 'botpres', zbotpres )
133
134      !                                         ! Mean density anomalie, temperature and salinity
135      ztemp = 0._wp
136      zsal  = 0._wp
137      DO jk = 1, jpkm1
138         DO jj = 1, jpj
139            DO ji = 1, jpi
140               zztmp = area(ji,jj) * fse3t(ji,jj,jk)
141               ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
142               zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
143            END DO
144         END DO
145      END DO
146      IF( .NOT.lk_vvl ) THEN
147         ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
148         zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
149      ENDIF
150      IF( lk_mpp ) THEN 
151         CALL mpp_sum( ztemp )
152         CALL mpp_sum( zsal  )
153      END IF
154      !
155      zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
156      ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
157      zsal  = zsal  / zvol                            ! Salinity of liquid seawater
158      !
159      CALL iom_put( 'masstot', zmass )
160      CALL iom_put( 'temptot', ztemp )
161      CALL iom_put( 'saltot' , zsal  )
162      !
163      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres )
164      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
165      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
166      !
167      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
168      !
169   END SUBROUTINE dia_ar5
170
171
172   SUBROUTINE dia_ar5_init
173      !!----------------------------------------------------------------------
174      !!                  ***  ROUTINE dia_ar5_init  ***
175      !!                   
176      !! ** Purpose :   initialization for AR5 diagnostic computation
177      !!----------------------------------------------------------------------
178      INTEGER  ::   inum
179      INTEGER  ::   ik
180      INTEGER  ::   ji, jj, jk  ! dummy loop indices
181      REAL(wp) ::   zztmp 
182      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
183      !!----------------------------------------------------------------------
184      !
185      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
186      !
187      CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
188      !                                      ! allocate dia_ar5 arrays
189      IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
190
191      area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:)
192
193      area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
194
195      vol0        = 0._wp
196      thick0(:,:) = 0._wp
197      DO jk = 1, jpkm1
198         vol0        = vol0        + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
199         thick0(:,:) = thick0(:,:) +    tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
200      END DO
201      IF( lk_mpp )   CALL mpp_sum( vol0 )
202     
203      CALL iom_open ( 'data_1m_salinity_nomask', inum )
204      CALL iom_get  ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,1), 1  )
205      CALL iom_get  ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,2), 12 )
206      CALL iom_close( inum )
207      sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
208      sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
209      IF( ln_zps ) THEN               ! z-coord. partial steps
210         DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
211            DO ji = 1, jpi
212               ik = mbkt(ji,jj)
213               IF( ik > 1 ) THEN
214                  zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
215                  sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
216               ENDIF
217            END DO
218         END DO
219      ENDIF
220      !
221      CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
222      !
223      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
224      !
225   END SUBROUTINE dia_ar5_init
226
227#else
228   !!----------------------------------------------------------------------
229   !!   Default option :                                         NO diaar5
230   !!----------------------------------------------------------------------
231   LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE.   ! coupled flag
232CONTAINS
233   SUBROUTINE dia_ar5_init    ! Dummy routine
234   END SUBROUTINE dia_ar5_init
235   SUBROUTINE dia_ar5( kt )   ! Empty routine
236      INTEGER ::   kt
237      WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
238   END SUBROUTINE dia_ar5
239#endif
240
241   !!======================================================================
242END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.