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 @ 12911

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

debug: gdept has to be defined over the whole water column

  • 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 (zgdept: needed to use the substitute)
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         DO jk = 1, jpk
135            zgdept(:,:,jk) = gdept(:,:,jk,Kmm)
136         END DO
137         CALL eos( ztsn, zrhd, zgdept)                       ! now in situ density using initial salinity
138         !
139         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
140         DO jk = 1, jpkm1
141            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * zrhd(:,:,jk)
142         END DO
143         IF( ln_linssh ) THEN
144            IF( ln_isfcav ) THEN
145               DO ji = 1, jpi
146                  DO jj = 1, jpj
147                     iks = mikt(ji,jj)
148                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * zrhd(ji,jj,iks) + riceload(ji,jj)
149                  END DO
150               END DO
151            ELSE
152               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
153            END IF
154!!gm
155!!gm   riceload should be added in both ln_linssh=T or F, no?
156!!gm
157         END IF
158         !                                         
159         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
160         zssh_steric = - zarho / area_tot
161         CALL iom_put( 'sshthster', zssh_steric )
162     
163         !                                         ! steric sea surface height
164         CALL eos( ts(:,:,:,:,Kmm), zrhd, zrhop, zgdept )                 ! now in situ and potential density
165         zrhop(:,:,jpk) = 0._wp
166         CALL iom_put( 'rhop', zrhop )
167         !
168         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
169         DO jk = 1, jpkm1
170            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * zrhd(:,:,jk)
171         END DO
172         IF( ln_linssh ) THEN
173            IF ( ln_isfcav ) THEN
174               DO ji = 1,jpi
175                  DO jj = 1,jpj
176                     iks = mikt(ji,jj)
177                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * zrhd(ji,jj,iks) + riceload(ji,jj)
178                  END DO
179               END DO
180            ELSE
181               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
182            END IF
183         END IF
184         !   
185         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
186         zssh_steric = - zarho / area_tot
187         CALL iom_put( 'sshsteric', zssh_steric )
188         !                                         ! ocean bottom pressure
189         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
190         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) )
191         CALL iom_put( 'botpres', zbotpres )
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_11_11( 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( ts(:,:,jk,jp_tem,Kmm), 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_11_11( 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_11_11( 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_11_11( 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( zrhd      , zrhop    )
303        DEALLOCATE( ztsn                 )
304      ENDIF
305      !
306      IF( ln_timing )   CALL timing_stop('dia_ar5')
307      !
308   END SUBROUTINE dia_ar5
309
310
311   SUBROUTINE dia_ar5_hst( ktra, cptr, puflx, pvflx ) 
312      !!----------------------------------------------------------------------
313      !!                    ***  ROUTINE dia_ar5_htr ***
314      !!----------------------------------------------------------------------
315      !! Wrapper for heat transport calculations
316      !! Called from all advection and/or diffusion routines
317      !!----------------------------------------------------------------------
318      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
319      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
320      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: puflx  ! u-flux of advection/diffusion
321      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pvflx  ! v-flux of advection/diffusion
322      !
323      INTEGER    ::  ji, jj, jk
324      REAL(wp), DIMENSION(jpi,jpj)  :: z2d
325
326   
327      z2d(:,:) = puflx(:,:,1) 
328      DO_3D_00_00( 1, jpkm1 )
329         z2d(ji,jj) = z2d(ji,jj) + puflx(ji,jj,jk) 
330      END_3D
331       CALL lbc_lnk( 'diaar5', z2d, 'U', -1. )
332       IF( cptr == 'adv' ) THEN
333          IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in i-direction
334          IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * z2d )  ! advective salt transport in i-direction
335       ENDIF
336       IF( cptr == 'ldf' ) THEN
337          IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in i-direction
338          IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * z2d ) ! diffusive salt transport in i-direction
339       ENDIF
340       !
341       z2d(:,:) = pvflx(:,:,1) 
342       DO_3D_00_00( 1, jpkm1 )
343          z2d(ji,jj) = z2d(ji,jj) + pvflx(ji,jj,jk) 
344       END_3D
345       CALL lbc_lnk( 'diaar5', z2d, 'V', -1. )
346       IF( cptr == 'adv' ) THEN
347          IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in j-direction
348          IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * z2d )  ! advective salt transport in j-direction
349       ENDIF
350       IF( cptr == 'ldf' ) THEN
351          IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in j-direction
352          IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * z2d ) ! diffusive salt transport in j-direction
353       ENDIF
354         
355   END SUBROUTINE dia_ar5_hst
356
357
358   SUBROUTINE dia_ar5_init
359      !!----------------------------------------------------------------------
360      !!                  ***  ROUTINE dia_ar5_init  ***
361      !!                   
362      !! ** Purpose :   initialization for AR5 diagnostic computation
363      !!----------------------------------------------------------------------
364      INTEGER  ::   inum
365      INTEGER  ::   ik, idep
366      INTEGER  ::   ji, jj, jk  ! dummy loop indices
367      REAL(wp) ::   zztmp 
368      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
369      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
370      !
371      !!----------------------------------------------------------------------
372      !
373      l_ar5 = .FALSE.
374      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
375         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
376         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) L_ar5 = .TRUE.
377 
378      IF( l_ar5 ) THEN
379         !
380         !                                      ! allocate dia_ar5 arrays
381         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
382
383         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
384
385         ALLOCATE( zvol0(jpi,jpj) )
386         zvol0 (:,:) = 0._wp
387         thick0(:,:) = 0._wp
388         DO_3D_11_11( 1, jpkm1 )
389            idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
390            zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
391            thick0(ji,jj) = thick0(ji,jj) +  idep   
392         END_3D
393         vol0 = glob_sum( 'diaar5', zvol0 )
394         DEALLOCATE( zvol0 )
395
396         IF( iom_use( 'sshthster' ) ) THEN
397            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
398            CALL iom_open ( 'sali_ref_clim_monthly', inum )
399            CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
400            CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
401            CALL iom_close( inum )
402
403            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
404            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
405            IF( ln_zps ) THEN               ! z-coord. partial steps
406               DO_2D_11_11
407                  ik = mbkt(ji,jj)
408                  IF( ik > 1 ) THEN
409                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
410                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
411                  ENDIF
412               END_2D
413            ENDIF
414            !
415            DEALLOCATE( zsaldta )
416         ENDIF
417         !
418      ENDIF
419      !
420   END SUBROUTINE dia_ar5_init
421
422   !!======================================================================
423END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.