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

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIA/diaar5.F90 @ 12193

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

2019/dev_r11943_MERGE_2019: Merge in dev_r12072_TOP-01_ENHANCE-11_cethe

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