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 NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DIA/diaar5.F90 @ 14219

Last change on this file since 14219 was 14219, checked in by mcastril, 4 years ago

Add Mixed Precision support by Oriol Tintó

  • Property svn:keywords set to Id
File size: 19.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   !!   dia_ar5       : AR5 diagnostics
10   !!   dia_ar5_init  : initialisation of AR5 diagnostics
11   !!----------------------------------------------------------------------
12   USE oce            ! ocean dynamics and active tracers
13   USE dom_oce        ! ocean space and time domain
14   USE eosbn2         ! equation of state                (eos_bn2 routine)
15   USE phycst         ! physical constant
16   USE in_out_manager  ! I/O manager
17   USE zdfddm
18   USE zdf_oce
19   !
20   USE lib_mpp        ! distribued memory computing library
21   USE iom            ! I/O manager library
22   USE fldread        ! type FLD_N
23   USE timing         ! preformance summary
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   dia_ar5        ! routine called in step.F90 module
29   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
30   PUBLIC   dia_ar5_hst    ! heat/salt transport
31
32   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
33   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
34   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:  ) ::   thick0       ! ocean thickness (interior domain)
35   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity
36   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   hstr_adv, hstr_ldf
37
38   LOGICAL  :: l_ar5
39
40   !! * Substitutions
41#  include "do_loop_substitute.h90"
42#  include "domzgr_substitute.h90"
43#  include "single_precision_substitute.h90"
44   !!----------------------------------------------------------------------
45   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
46   !! $Id$
47   !! Software governed by the CeCILL license (see ./LICENSE)
48   !!----------------------------------------------------------------------
49CONTAINS
50
51   FUNCTION dia_ar5_alloc()
52      !!----------------------------------------------------------------------
53      !!                    ***  ROUTINE dia_ar5_alloc  ***
54      !!----------------------------------------------------------------------
55      INTEGER :: dia_ar5_alloc
56      !!----------------------------------------------------------------------
57      !
58      ALLOCATE( thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , &
59         &      hstr_adv(jpi,jpj,jpts,2), hstr_ldf(jpi,jpj,jpts,2), STAT=dia_ar5_alloc )
60      !
61      CALL mpp_sum ( 'diaar5', dia_ar5_alloc )
62      IF( dia_ar5_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate arrays' )
63      !
64   END FUNCTION dia_ar5_alloc
65
66
67   SUBROUTINE dia_ar5( kt, Kmm )
68      !!----------------------------------------------------------------------
69      !!                    ***  ROUTINE dia_ar5  ***
70      !!
71      !! ** Purpose :   compute and output some AR5 diagnostics
72      !!----------------------------------------------------------------------
73      !
74      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
75      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
76      !
77      INTEGER  ::   ji, jj, jk, iks, ikb                      ! dummy loop arguments
78      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass, zsst
79      REAL(wp) ::   zaw, zbw, zrw
80      !
81      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
82      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z2d, zpe                   ! 2D workspace
83      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: z3d, zrhd, ztpot, zgdept   ! 3D workspace (zgdept: needed to use the substitute)
84      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
85
86      !!--------------------------------------------------------------------
87      IF( ln_timing )   CALL timing_start('dia_ar5')
88
89      IF( kt == nit000 )     CALL dia_ar5_init
90
91      IF( l_ar5 ) THEN
92         ALLOCATE( zarea_ssh(jpi,jpj), zbotpres(jpi,jpj), z2d(jpi,jpj) )
93         ALLOCATE( zrhd(jpi,jpj,jpk) )
94         ALLOCATE( ztsn(jpi,jpj,jpk,jpts) )
95         zarea_ssh(:,:) = e1e2t(:,:) * ssh(:,:,Kmm)
96      ENDIF
97      !
98      CALL iom_put( 'e2u'      , e2u  (:,:) )
99      CALL iom_put( 'e1v'      , e1v  (:,:) )
100      CALL iom_put( 'areacello', e1e2t(:,:) )
101      !
102      IF( iom_use( 'volcello' ) .OR. iom_use( 'masscello' )  ) THEN
103         zrhd(:,:,jpk) = 0._wp        ! ocean volume ; rhd is used as workspace
104         DO jk = 1, jpkm1
105            zrhd(:,:,jk) = e1e2t(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
106         END DO
107         DO jk = 1, jpk
108            z3d(:,:,jk) =  rho0 * e3t(:,:,jk,Kmm) * tmask(:,:,jk)
109         END DO
110         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
111         CALL iom_put( 'masscello' , z3d (:,:,:) )   ! ocean mass
112      ENDIF
113      !
114      IF( iom_use( 'e3tb' ) )  THEN    ! bottom layer thickness
115         DO_2D( 1, 1, 1, 1 )
116            ikb = mbkt(ji,jj)
117            z2d(ji,jj) = e3t(ji,jj,ikb,Kmm)
118         END_2D
119         CALL iom_put( 'e3tb', z2d )
120      ENDIF
121      !
122      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN
123         !                                         ! total volume of liquid seawater
124         zvolssh = glob_sum( 'diaar5', zarea_ssh(:,:) )
125         zvol    = vol0 + zvolssh
126
127         CALL iom_put( 'voltot', zvol               )
128         CALL iom_put( 'sshtot', zvolssh / area_tot )
129         CALL iom_put( 'sshdyn', ssh(:,:,Kmm) - (zvolssh / area_tot) )
130         !
131      ENDIF
132
133      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN
134         !
135         ztsn(:,:,:,jp_tem) = ts(:,:,:,jp_tem,Kmm)                    ! thermosteric ssh
136         ztsn(:,:,:,jp_sal) = sn0(:,:,:)
137         ALLOCATE( zgdept(jpi,jpj,jpk) )
138         DO jk = 1, jpk
139            zgdept(:,:,jk) = gdept(:,:,jk,Kmm)
140         END DO
141         CALL eos( ztsn, zrhd, zgdept)                       ! now in situ density using initial salinity
142         !
143         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
144         DO jk = 1, jpkm1
145            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * zrhd(:,:,jk)
146         END DO
147         IF( ln_linssh ) THEN
148            IF( ln_isfcav ) THEN
149               DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
150                  iks = mikt(ji,jj)
151                  zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * zrhd(ji,jj,iks) + riceload(ji,jj)
152               END_2D
153            ELSE
154               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
155            END IF
156!!gm
157!!gm   riceload should be added in both ln_linssh=T or F, no?
158!!gm
159         END IF
160         !
161         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) )
162         zssh_steric = - zarho / area_tot
163         CALL iom_put( 'sshthster', zssh_steric )
164
165         !                                         ! steric sea surface height
166         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
167         DO jk = 1, jpkm1
168            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * rhd(:,:,jk)
169         END DO
170         IF( ln_linssh ) THEN
171            IF ( ln_isfcav ) THEN
172               DO ji = 1,jpi
173                  DO jj = 1,jpj
174                     iks = mikt(ji,jj)
175                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * rhd(ji,jj,iks) + riceload(ji,jj)
176                  END DO
177               END DO
178            ELSE
179               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * rhd(:,:,1)
180            END IF
181         END IF
182         !
183         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) )
184         zssh_steric = - zarho / area_tot
185         CALL iom_put( 'sshsteric', zssh_steric )
186         !                                         ! ocean bottom pressure
187         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
188         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) )
189         CALL iom_put( 'botpres', zbotpres )
190         !
191         DEALLOCATE( zgdept )
192         !
193      ENDIF
194
195      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN
196          !                                         ! Mean density anomalie, temperature and salinity
197          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
198          DO_3D( 1, 1, 1, 1, 1, jpkm1 )
199             zztmp = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm)
200             ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * ts(ji,jj,jk,jp_tem,Kmm)
201             ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * ts(ji,jj,jk,jp_sal,Kmm)
202          END_3D
203
204          IF( ln_linssh ) THEN
205            IF( ln_isfcav ) THEN
206               DO ji = 1, jpi
207                  DO jj = 1, jpj
208                     iks = mikt(ji,jj)
209                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_tem,Kmm)
210                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_sal,Kmm)
211                  END DO
212               END DO
213            ELSE
214               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * ts(:,:,1,jp_tem,Kmm)
215               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * ts(:,:,1,jp_sal,Kmm)
216            END IF
217         ENDIF
218         !
219         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) )
220         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) )
221         zmass = rho0 * ( zarho + zvol )
222         !
223         CALL iom_put( 'masstot', zmass )
224         CALL iom_put( 'temptot', ztemp / zvol )
225         CALL iom_put( 'saltot' , zsal  / zvol )
226         !
227      ENDIF
228
229      IF( ln_teos10 ) THEN        ! ! potential temperature (TEOS-10 case)
230         IF( iom_use( 'toce_pot') .OR. iom_use( 'temptot_pot' ) .OR. iom_use( 'sst_pot' )  &
231                                  .OR. iom_use( 'ssttot' ) .OR.  iom_use( 'tosmint_pot' ) ) THEN
232            !
233            ALLOCATE( ztpot(jpi,jpj,jpk) )
234            ztpot(:,:,jpk) = 0._wp
235            DO jk = 1, jpkm1
236               ztpot(:,:,jk) = eos_pt_from_ct( CASTWP(ts(:,:,jk,jp_tem,Kmm)), CASTWP(ts(:,:,jk,jp_sal,Kmm)) )
237            END DO
238            !
239            CALL iom_put( 'toce_pot', ztpot(:,:,:) )  ! potential temperature (TEOS-10 case)
240            CALL iom_put( 'sst_pot' , ztpot(:,:,1) )  ! surface temperature
241            !
242            IF( iom_use( 'temptot_pot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
243               z2d(:,:) = 0._wp
244               DO jk = 1, jpkm1
245                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t(:,:,jk,Kmm) * ztpot(:,:,jk)
246               END DO
247               ztemp = glob_sum( 'diaar5', z2d(:,:)  )
248               CALL iom_put( 'temptot_pot', ztemp / zvol )
249             ENDIF
250             !
251             IF( iom_use( 'ssttot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
252               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  )
253               CALL iom_put( 'ssttot', zsst / area_tot )
254             ENDIF
255             ! Vertical integral of temperature
256             IF( iom_use( 'tosmint_pot') ) THEN
257               z2d(:,:) = 0._wp
258               DO_3D( 1, 1, 1, 1, 1, jpkm1 )
259                  z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk)
260               END_3D
261               CALL iom_put( 'tosmint_pot', z2d )
262            ENDIF
263            DEALLOCATE( ztpot )
264        ENDIF
265      ELSE
266         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
267            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * ts(:,:,1,jp_tem,Kmm) )
268            CALL iom_put('ssttot', zsst / area_tot )
269         ENDIF
270      ENDIF
271
272      IF( iom_use( 'tnpeo' )) THEN   
273        ! Work done against stratification by vertical mixing
274        ! Exclude points where rn2 is negative as convection kicks in here and
275        ! work is not being done against stratification
276         ALLOCATE( zpe(jpi,jpj) )
277         zpe(:,:) = 0._wp
278         IF( ln_zdfddm ) THEN
279            DO_3D( 1, 1, 1, 1, 2, jpk )
280               IF( rn2(ji,jj,jk) > 0._wp ) THEN
281                  zrw = ( gdept(ji,jj,jk,Kmm) - gdepw(ji,jj,jk,Kmm) ) / e3w(ji,jj,jk,Kmm)
282                  !
283                  zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
284                  zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
285                  !
286                  zpe(ji, jj) = zpe(ji,jj)   &
287                     &        -  grav * (  avt(ji,jj,jk) * zaw * (ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) )  &
288                     &                   - avs(ji,jj,jk) * zbw * (ts(ji,jj,jk-1,jp_sal,Kmm) - ts(ji,jj,jk,jp_sal,Kmm) ) )
289               ENDIF
290            END_3D
291          ELSE
292            DO_3D( 1, 1, 1, 1, 1, jpk )
293               zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rho0 * e3w(ji,jj,jk,Kmm)
294            END_3D
295         ENDIF
296          CALL iom_put( 'tnpeo', zpe )
297          DEALLOCATE( zpe )
298      ENDIF
299
300      IF( l_ar5 ) THEN
301        DEALLOCATE( zarea_ssh , zbotpres, z2d )
302        DEALLOCATE( ztsn                 )
303      ENDIF
304      !
305      IF( ln_timing )   CALL timing_stop('dia_ar5')
306      !
307   END SUBROUTINE dia_ar5
308
309   ! TEMP: [tiling] These changes not necessary if using XIOS (subdomain support, will not output haloes)
310   SUBROUTINE dia_ar5_hst( ktra, cptr, puflx, pvflx )
311      !!----------------------------------------------------------------------
312      !!                    ***  ROUTINE dia_ar5_htr ***
313      !!----------------------------------------------------------------------
314      !! Wrapper for heat transport calculations
315      !! Called from all advection and/or diffusion routines
316      !!----------------------------------------------------------------------
317      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
318      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
319      REAL(wp), DIMENSION(A2D(nn_hls),jpk)    , INTENT(in)   :: puflx  ! u-flux of advection/diffusion
320      REAL(wp), DIMENSION(A2D(nn_hls),jpk)    , INTENT(in)   :: pvflx  ! v-flux of advection/diffusion
321      !
322      INTEGER    ::  ji, jj, jk
323
324      IF( cptr /= 'adv' .AND. cptr /= 'ldf' ) RETURN
325      IF( ktra /= jp_tem .AND. ktra /= jp_sal ) RETURN
326
327      IF( cptr == 'adv' ) THEN
328         DO_2D( 0, 0, 0, 0 )
329            hstr_adv(ji,jj,ktra,1) = puflx(ji,jj,1)
330            hstr_adv(ji,jj,ktra,2) = pvflx(ji,jj,1)
331         END_2D
332         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
333            hstr_adv(ji,jj,ktra,1) = hstr_adv(ji,jj,ktra,1) + puflx(ji,jj,jk)
334            hstr_adv(ji,jj,ktra,2) = hstr_adv(ji,jj,ktra,2) + pvflx(ji,jj,jk)
335         END_3D
336      ELSE IF( cptr == 'ldf' ) THEN
337         DO_2D( 0, 0, 0, 0 )
338            hstr_ldf(ji,jj,ktra,1) = puflx(ji,jj,1)
339            hstr_ldf(ji,jj,ktra,2) = pvflx(ji,jj,1)
340         END_2D
341         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
342            hstr_ldf(ji,jj,ktra,1) = hstr_ldf(ji,jj,ktra,1) + puflx(ji,jj,jk)
343            hstr_ldf(ji,jj,ktra,2) = hstr_ldf(ji,jj,ktra,2) + pvflx(ji,jj,jk)
344         END_3D
345      ENDIF
346
347      IF( ntile == 0 .OR. ntile == nijtile ) THEN
348         IF( cptr == 'adv' ) THEN
349            IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,1) )  ! advective heat transport in i-direction
350            IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * hstr_adv(:,:,ktra,1) )  ! advective salt transport in i-direction
351            IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * hstr_adv(:,:,ktra,2) )  ! advective heat transport in j-direction
352            IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * hstr_adv(:,:,ktra,2) )  ! advective salt transport in j-direction
353         ENDIF
354         IF( cptr == 'ldf' ) THEN
355            IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,1) ) ! diffusive heat transport in i-direction
356            IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * hstr_ldf(:,:,ktra,1) ) ! diffusive salt transport in i-direction
357            IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * hstr_ldf(:,:,ktra,2) ) ! diffusive heat transport in j-direction
358            IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * hstr_ldf(:,:,ktra,2) ) ! diffusive salt transport in j-direction
359         ENDIF
360      ENDIF
361
362   END SUBROUTINE dia_ar5_hst
363
364
365   SUBROUTINE dia_ar5_init
366      !!----------------------------------------------------------------------
367      !!                  ***  ROUTINE dia_ar5_init  ***
368      !!
369      !! ** Purpose :   initialization for AR5 diagnostic computation
370      !!----------------------------------------------------------------------
371      INTEGER  ::   inum
372      INTEGER  ::   ik, idep
373      INTEGER  ::   ji, jj, jk  ! dummy loop indices
374      REAL(wp) ::   zztmp
375      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
376      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0
377      !
378      !!----------------------------------------------------------------------
379      !
380      l_ar5 = .FALSE.
381      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  &
382         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &
383         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' ) .OR. &
384         &  iom_use( 'uadv_heattr' ) .OR. iom_use( 'udiff_heattr' ) .OR. &
385         &  iom_use( 'uadv_salttr' ) .OR. iom_use( 'udiff_salttr' ) .OR. &
386         &  iom_use( 'vadv_heattr' ) .OR. iom_use( 'vdiff_heattr' ) .OR. &
387         &  iom_use( 'vadv_salttr' ) .OR. iom_use( 'vdiff_salttr' ) .OR. &
388         &  iom_use( 'rhop' )  ) L_ar5 = .TRUE.
389
390      IF( l_ar5 ) THEN
391         !
392         !                                      ! allocate dia_ar5 arrays
393         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
394
395         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
396
397         ALLOCATE( zvol0(jpi,jpj) )
398         zvol0 (:,:) = 0._wp
399         thick0(:,:) = 0._wp
400         DO_3D( 1, 1, 1, 1, 1, jpkm1 )   ! interpolation of salinity at the last ocean level (i.e. the partial step)
401            idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
402            zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
403            thick0(ji,jj) = thick0(ji,jj) +  idep
404         END_3D
405         vol0 = glob_sum( 'diaar5', zvol0 )
406         DEALLOCATE( zvol0 )
407
408         IF( iom_use( 'sshthster' ) ) THEN
409            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
410            CALL iom_open ( 'sali_ref_clim_monthly', inum )
411            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,1), 1  )
412            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,2), 12 )
413            CALL iom_close( inum )
414
415            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )
416            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
417            IF( ln_zps ) THEN               ! z-coord. partial steps
418               DO_2D( 1, 1, 1, 1 )          ! interpolation of salinity at the last ocean level (i.e. the partial step)
419                  ik = mbkt(ji,jj)
420                  IF( ik > 1 ) THEN
421                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
422                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
423                  ENDIF
424               END_2D
425            ENDIF
426            !
427            DEALLOCATE( zsaldta )
428         ENDIF
429         !
430      ENDIF
431      !
432   END SUBROUTINE dia_ar5_init
433
434   !!======================================================================
435END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.