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

Last change on this file since 12712 was 12712, checked in by davestorkey, 8 months ago

branches/2020/r4.0-HEAD_ticket2425 : first guess solution

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