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

source: NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser/src/OCE/DIA/diaar5.F90 @ 12928

Last change on this file since 12928 was 12928, checked in by smueller, 4 years ago

Synchronizing with /NEMO/trunk@12925 (ticket #2170)

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