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/trunk/src/OCE/DIA – NEMO

source: NEMO/trunk/src/OCE/DIA/diaar5.F90 @ 13497

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

re-introduce comments that have been erased by loops transformation see #2525

  • Property svn:keywords set to Id
File size: 18.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
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, 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) )
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( 1, 1, 1, 1 )
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         ALLOCATE( zgdept(jpi,jpj,jpk) )
135         DO jk = 1, jpk
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_2D( nn_hls, nn_hls, nn_hls, nn_hls )
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_2D
150            ELSE
151               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * zrhd(:,:,1)
152            END IF
153!!gm
154!!gm   riceload should be added in both ln_linssh=T or F, no?
155!!gm
156         END IF
157         !                                         
158         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
159         zssh_steric = - zarho / area_tot
160         CALL iom_put( 'sshthster', zssh_steric )
161     
162         !                                         ! steric sea surface height
163         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
164         DO jk = 1, jpkm1
165            zbotpres(:,:) = zbotpres(:,:) + e3t(:,:,jk,Kmm) * rhd(:,:,jk)
166         END DO
167         IF( ln_linssh ) THEN
168            IF ( ln_isfcav ) THEN
169               DO ji = 1,jpi
170                  DO jj = 1,jpj
171                     iks = mikt(ji,jj)
172                     zbotpres(ji,jj) = zbotpres(ji,jj) + ssh(ji,jj,Kmm) * rhd(ji,jj,iks) + riceload(ji,jj)
173                  END DO
174               END DO
175            ELSE
176               zbotpres(:,:) = zbotpres(:,:) + ssh(:,:,Kmm) * rhd(:,:,1)
177            END IF
178         END IF
179         !   
180         zarho = glob_sum( 'diaar5', e1e2t(:,:) * zbotpres(:,:) ) 
181         zssh_steric = - zarho / area_tot
182         CALL iom_put( 'sshsteric', zssh_steric )
183         !                                         ! ocean bottom pressure
184         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
185         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) )
186         CALL iom_put( 'botpres', zbotpres )
187         !
188         DEALLOCATE( zgdept )
189         !
190      ENDIF
191
192      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
193          !                                         ! Mean density anomalie, temperature and salinity
194          ztsn(:,:,:,:) = 0._wp                    ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
195          DO_3D( 1, 1, 1, 1, 1, jpkm1 )
196             zztmp = e1e2t(ji,jj) * e3t(ji,jj,jk,Kmm)
197             ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zztmp * ts(ji,jj,jk,jp_tem,Kmm)
198             ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zztmp * ts(ji,jj,jk,jp_sal,Kmm)
199          END_3D
200
201          IF( ln_linssh ) THEN
202            IF( ln_isfcav ) THEN
203               DO ji = 1, jpi
204                  DO jj = 1, jpj
205                     iks = mikt(ji,jj)
206                     ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_tem,Kmm) 
207                     ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * ts(ji,jj,iks,jp_sal,Kmm) 
208                  END DO
209               END DO
210            ELSE
211               ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * ts(:,:,1,jp_tem,Kmm) 
212               ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * ts(:,:,1,jp_sal,Kmm) 
213            END IF
214         ENDIF
215         !
216         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) )
217         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) )
218         zmass = rho0 * ( zarho + zvol )     
219         !
220         CALL iom_put( 'masstot', zmass )
221         CALL iom_put( 'temptot', ztemp / zvol )
222         CALL iom_put( 'saltot' , zsal  / zvol )
223         !
224      ENDIF     
225
226      IF( ln_teos10 ) THEN        ! ! potential temperature (TEOS-10 case)
227         IF( iom_use( 'toce_pot') .OR. iom_use( 'temptot_pot' ) .OR. iom_use( 'sst_pot' )  &
228                                  .OR. iom_use( 'ssttot' ) .OR.  iom_use( 'tosmint_pot' ) ) THEN
229            !
230            ALLOCATE( ztpot(jpi,jpj,jpk) )
231            ztpot(:,:,jpk) = 0._wp
232            DO jk = 1, jpkm1
233               ztpot(:,:,jk) = eos_pt_from_ct( ts(:,:,jk,jp_tem,Kmm), ts(:,:,jk,jp_sal,Kmm) )
234            END DO
235            !
236            CALL iom_put( 'toce_pot', ztpot(:,:,:) )  ! potential temperature (TEOS-10 case)
237            CALL iom_put( 'sst_pot' , ztpot(:,:,1) )  ! surface temperature
238            !
239            IF( iom_use( 'temptot_pot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
240               z2d(:,:) = 0._wp
241               DO jk = 1, jpkm1
242                 z2d(:,:) = z2d(:,:) + e1e2t(:,:) * e3t(:,:,jk,Kmm) * ztpot(:,:,jk)
243               END DO
244               ztemp = glob_sum( 'diaar5', z2d(:,:)  ) 
245               CALL iom_put( 'temptot_pot', ztemp / zvol )
246             ENDIF
247             !
248             IF( iom_use( 'ssttot' ) ) THEN   ! Output potential temperature in case we use TEOS-10
249               zsst = glob_sum( 'diaar5',  e1e2t(:,:) * ztpot(:,:,1)  ) 
250               CALL iom_put( 'ssttot', zsst / area_tot )
251             ENDIF
252             ! Vertical integral of temperature
253             IF( iom_use( 'tosmint_pot') ) THEN
254               z2d(:,:) = 0._wp
255               DO_3D( 1, 1, 1, 1, 1, jpkm1 )
256                  z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk)
257               END_3D
258               CALL iom_put( 'tosmint_pot', z2d ) 
259            ENDIF
260            DEALLOCATE( ztpot )
261        ENDIF
262      ELSE       
263         IF( iom_use('ssttot') ) THEN   ! Output sst in case we use EOS-80
264            zsst  = glob_sum( 'diaar5', e1e2t(:,:) * ts(:,:,1,jp_tem,Kmm) )
265            CALL iom_put('ssttot', zsst / area_tot )
266         ENDIF
267      ENDIF
268
269      IF( iom_use( 'tnpeo' )) THEN   
270        ! Work done against stratification by vertical mixing
271        ! Exclude points where rn2 is negative as convection kicks in here and
272        ! work is not being done against stratification
273         ALLOCATE( zpe(jpi,jpj) )
274         zpe(:,:) = 0._wp
275         IF( ln_zdfddm ) THEN
276            DO_3D( 1, 1, 1, 1, 2, jpk )
277               IF( rn2(ji,jj,jk) > 0._wp ) THEN
278                  zrw = ( gdept(ji,jj,jk,Kmm) - gdepw(ji,jj,jk,Kmm) ) / e3w(ji,jj,jk,Kmm)
279                  !
280                  zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
281                  zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
282                  !
283                  zpe(ji, jj) = zpe(ji,jj)   &
284                     &        -  grav * (  avt(ji,jj,jk) * zaw * (ts(ji,jj,jk-1,jp_tem,Kmm) - ts(ji,jj,jk,jp_tem,Kmm) )  &
285                     &                   - avs(ji,jj,jk) * zbw * (ts(ji,jj,jk-1,jp_sal,Kmm) - ts(ji,jj,jk,jp_sal,Kmm) ) )
286               ENDIF
287            END_3D
288          ELSE
289            DO_3D( 1, 1, 1, 1, 1, jpk )
290               zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rho0 * e3w(ji,jj,jk,Kmm)
291            END_3D
292         ENDIF
293          CALL iom_put( 'tnpeo', zpe )
294          DEALLOCATE( zpe )
295      ENDIF
296
297      IF( l_ar5 ) THEN
298        DEALLOCATE( zarea_ssh , zbotpres, z2d )
299        DEALLOCATE( ztsn                 )
300      ENDIF
301      !
302      IF( ln_timing )   CALL timing_stop('dia_ar5')
303      !
304   END SUBROUTINE dia_ar5
305
306
307   SUBROUTINE dia_ar5_hst( ktra, cptr, puflx, pvflx ) 
308      !!----------------------------------------------------------------------
309      !!                    ***  ROUTINE dia_ar5_htr ***
310      !!----------------------------------------------------------------------
311      !! Wrapper for heat transport calculations
312      !! Called from all advection and/or diffusion routines
313      !!----------------------------------------------------------------------
314      INTEGER                         , INTENT(in )  :: ktra  ! tracer index
315      CHARACTER(len=3)                , INTENT(in)   :: cptr  ! transport type  'adv'/'ldf'
316      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: puflx  ! u-flux of advection/diffusion
317      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pvflx  ! v-flux of advection/diffusion
318      !
319      INTEGER    ::  ji, jj, jk
320      REAL(wp), DIMENSION(jpi,jpj)  :: z2d
321
322   
323      z2d(:,:) = puflx(:,:,1) 
324      DO_3D( 0, 0, 0, 0, 1, jpkm1 )
325         z2d(ji,jj) = z2d(ji,jj) + puflx(ji,jj,jk) 
326      END_3D
327       CALL lbc_lnk( 'diaar5', z2d, 'U', -1.0_wp )
328       IF( cptr == 'adv' ) THEN
329          IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in i-direction
330          IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * z2d )  ! advective salt transport in i-direction
331       ENDIF
332       IF( cptr == 'ldf' ) THEN
333          IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in i-direction
334          IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * z2d ) ! diffusive salt transport in i-direction
335       ENDIF
336       !
337       z2d(:,:) = pvflx(:,:,1) 
338       DO_3D( 0, 0, 0, 0, 1, jpkm1 )
339          z2d(ji,jj) = z2d(ji,jj) + pvflx(ji,jj,jk) 
340       END_3D
341       CALL lbc_lnk( 'diaar5', z2d, 'V', -1.0_wp )
342       IF( cptr == 'adv' ) THEN
343          IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in j-direction
344          IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * z2d )  ! advective salt transport in j-direction
345       ENDIF
346       IF( cptr == 'ldf' ) THEN
347          IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in j-direction
348          IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * z2d ) ! diffusive salt transport in j-direction
349       ENDIF
350         
351   END SUBROUTINE dia_ar5_hst
352
353
354   SUBROUTINE dia_ar5_init
355      !!----------------------------------------------------------------------
356      !!                  ***  ROUTINE dia_ar5_init  ***
357      !!                   
358      !! ** Purpose :   initialization for AR5 diagnostic computation
359      !!----------------------------------------------------------------------
360      INTEGER  ::   inum
361      INTEGER  ::   ik, idep
362      INTEGER  ::   ji, jj, jk  ! dummy loop indices
363      REAL(wp) ::   zztmp 
364      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
365      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     ::   zvol0     
366      !
367      !!----------------------------------------------------------------------
368      !
369      l_ar5 = .FALSE.
370      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
371         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
372         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' ) .OR. &
373         &  iom_use( 'rhop' )  ) L_ar5 = .TRUE.
374 
375      IF( l_ar5 ) THEN
376         !
377         !                                      ! allocate dia_ar5 arrays
378         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
379
380         area_tot  = glob_sum( 'diaar5', e1e2t(:,:) )
381
382         ALLOCATE( zvol0(jpi,jpj) )
383         zvol0 (:,:) = 0._wp
384         thick0(:,:) = 0._wp
385         DO_3D( 1, 1, 1, 1, 1, jpkm1 )   ! interpolation of salinity at the last ocean level (i.e. the partial step)
386            idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
387            zvol0 (ji,jj) = zvol0 (ji,jj) +  idep * e1e2t(ji,jj)
388            thick0(ji,jj) = thick0(ji,jj) +  idep   
389         END_3D
390         vol0 = glob_sum( 'diaar5', zvol0 )
391         DEALLOCATE( zvol0 )
392
393         IF( iom_use( 'sshthster' ) ) THEN
394            ALLOCATE( zsaldta(jpi,jpj,jpk,jpts) )
395            CALL iom_open ( 'sali_ref_clim_monthly', inum )
396            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,1), 1  )
397            CALL iom_get  ( inum, jpdom_global, 'vosaline' , zsaldta(:,:,:,2), 12 )
398            CALL iom_close( inum )
399
400            sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )       
401            sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
402            IF( ln_zps ) THEN               ! z-coord. partial steps
403               DO_2D( 1, 1, 1, 1 )          ! interpolation of salinity at the last ocean level (i.e. the partial step)
404                  ik = mbkt(ji,jj)
405                  IF( ik > 1 ) THEN
406                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
407                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
408                  ENDIF
409               END_2D
410            ENDIF
411            !
412            DEALLOCATE( zsaldta )
413         ENDIF
414         !
415      ENDIF
416      !
417   END SUBROUTINE dia_ar5_init
418
419   !!======================================================================
420END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.