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_r12377_KERNEL-06_techene_e3/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/DIA/diaar5.F90 @ 12731

Last change on this file since 12731 was 12731, checked in by techene, 4 years ago

replace h. and gde. in case key_qco is activated - quick and dirty

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