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/r4.0-HEAD_ticket2425/src/OCE/DIA – NEMO

source: NEMO/branches/2020/r4.0-HEAD_ticket2425/src/OCE/DIA/diaar5.F90 @ 12727

Last change on this file since 12727 was 12727, checked in by davestorkey, 4 years ago

branches/2020/r4.0-HEAD_ticket2425 : Second version.
Remove unnecessary second eos call from dia_ar5 and move iom_put for rhop to dia_wri.

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