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_r12953_ENHANCE-10_acc_fix_traqsr/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r12953_ENHANCE-10_acc_fix_traqsr/src/OCE/DIA/diaar5.F90 @ 13121

Last change on this file since 13121 was 13121, checked in by acc, 4 years ago

2020/dev_r12953_ENHANCE-10_acc_fix_traqsr. Merge in trunk changes from 12953 to current HEAD (13115). Fully SETTE tested

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