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 trunk/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90 @ 7698

Last change on this file since 7698 was 7698, checked in by mocavero, 7 years ago

update trunk with OpenMP parallelization

  • Property svn:keywords set to Id
File size: 19.6 KB
RevLine 
[1756]1MODULE diaar5
2   !!======================================================================
3   !!                       ***  MODULE  diaar5  ***
4   !! AR5 diagnostics
5   !!======================================================================
[2528]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
[1756]8   !!----------------------------------------------------------------------
[2528]9   !!   dia_ar5       : AR5 diagnostics
10   !!   dia_ar5_init  : initialisation of AR5 diagnostics
[1756]11   !!----------------------------------------------------------------------
12   USE oce            ! ocean dynamics and active tracers
13   USE dom_oce        ! ocean space and time domain
[2528]14   USE eosbn2         ! equation of state                (eos_bn2 routine)
[1756]15   USE lib_mpp        ! distribued memory computing library
16   USE iom            ! I/O manager library
[3294]17   USE timing         ! preformance summary
18   USE wrk_nemo       ! working arrays
[5253]19   USE fldread        ! type FLD_N
20   USE phycst         ! physical constant
21   USE in_out_manager  ! I/O manager
[7646]22   USE zdfddm
23   USE zdf_oce
[1756]24
25   IMPLICIT NONE
26   PRIVATE
27
[2528]28   PUBLIC   dia_ar5        ! routine called in step.F90 module
[2715]29   PUBLIC   dia_ar5_alloc  ! routine called in nemogcm.F90 module
[7646]30   PUBLIC   dia_ar5_hst    ! heat/salt transport
[1756]31
[2528]32   REAL(wp)                         ::   vol0         ! ocean volume (interior domain)
33   REAL(wp)                         ::   area_tot     ! total ocean surface (interior domain)
[2715]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
[7646]37
38   LOGICAL  :: l_ar5
[1756]39     
[7646]40   !! * Substitutions
41#  include "zdfddm_substitute.h90"
42#  include "vectopt_loop_substitute.h90"
[1756]43   !!----------------------------------------------------------------------
[2528]44   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
45   !! $Id$
46   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[1756]47   !!----------------------------------------------------------------------
48CONTAINS
49
[2715]50   FUNCTION dia_ar5_alloc()
51      !!----------------------------------------------------------------------
52      !!                    ***  ROUTINE dia_ar5_alloc  ***
53      !!----------------------------------------------------------------------
54      INTEGER :: dia_ar5_alloc
55      !!----------------------------------------------------------------------
56      !
57      ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
58      !
59      IF( lk_mpp             )   CALL mpp_sum ( dia_ar5_alloc )
60      IF( dia_ar5_alloc /= 0 )   CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
61      !
62   END FUNCTION dia_ar5_alloc
63
64
[1756]65   SUBROUTINE dia_ar5( kt )
66      !!----------------------------------------------------------------------
67      !!                    ***  ROUTINE dia_ar5  ***
68      !!
[2528]69      !! ** Purpose :   compute and output some AR5 diagnostics
[1756]70      !!----------------------------------------------------------------------
[2715]71      !
[1756]72      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
[2715]73      !
[1756]74      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments
75      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
[7646]76      REAL(wp) ::   zaw, zbw, zrw
[3294]77      !
78      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace
[7646]79      REAL(wp), POINTER, DIMENSION(:,:)     :: zpe                         ! 2D workspace
[3294]80      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace
81      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace
[1756]82      !!--------------------------------------------------------------------
[3294]83      IF( nn_timing == 1 )   CALL timing_start('dia_ar5')
84 
[7646]85      IF( kt == nit000 )     CALL dia_ar5_init
[1756]86
[7646]87      IF( l_ar5 ) THEN
88         CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres )
89         CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    )
90         CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 )
[7698]91!$OMP PARALLEL DO schedule(static) private(jj, ji)
92         DO jj = 1, jpj
93            DO ji = 1, jpi
94               zarea_ssh(ji,jj) = area(ji,jj) * sshn(ji,jj)
95            END DO
96         END DO
[7646]97      ENDIF
98      !
99      IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' )  .OR. iom_use( 'sshdyn' )  ) THEN   
100         !                                         ! total volume of liquid seawater
101         zvolssh = SUM( zarea_ssh(:,:) ) 
102         IF( lk_mpp )   CALL mpp_sum( zvolssh )
103         zvol = vol0 + zvolssh
[1756]104     
[7646]105         CALL iom_put( 'voltot', zvol               )
106         CALL iom_put( 'sshtot', zvolssh / area_tot )
107         CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
108         !
109      ENDIF
[1756]110
[7646]111      IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) THEN   
112         !                     
[7698]113!$OMP PARALLEL DO schedule(static) private(jk, jj, ji)
114         DO jk = 1, jpk
115            DO jj = 1, jpj
116               DO ji = 1, jpi
117                  ztsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem)                    ! thermosteric ssh
118                  ztsn(ji,jj,jk,jp_sal) = sn0(ji,jj,jk)
119               END DO
120            END DO
121         END DO
[7646]122         CALL eos( ztsn, zrhd, gdept_n(:,:,:) )                       ! now in situ density using initial salinity
123         !
[7698]124!$OMP PARALLEL
125!$OMP DO schedule(static) private(jj, ji)
126         DO jj = 1, jpj
127            DO ji = 1, jpi
128               zbotpres(ji,jj) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
129            END DO
130         END DO
[7646]131         DO jk = 1, jpkm1
[7698]132!$OMP DO schedule(static) private(jj, ji)
133            DO jj = 1, jpj
134               DO ji = 1, jpi
135                  zbotpres(ji,jj) = zbotpres(ji,jj) + e3t_n(ji,jj,jk) * zrhd(ji,jj,jk)
136               END DO
137            END DO
[7646]138         END DO
[7698]139!$OMP END PARALLEL
[7646]140         IF( ln_linssh ) THEN
141            IF( ln_isfcav ) THEN
[7698]142!$OMP PARALLEL DO schedule(static) private(jj, ji)
[7646]143               DO ji = 1, jpi
144                  DO jj = 1, jpj
145                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
146                  END DO
[5120]147               END DO
[7646]148            ELSE
[7698]149!$OMP PARALLEL DO schedule(static) private(jj, ji)
150               DO ji = 1, jpi
151                  DO jj = 1, jpj
152                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,1)
153                  END DO
154               END DO
[7646]155            END IF
[6140]156!!gm
157!!gm   riceload should be added in both ln_linssh=T or F, no?
158!!gm
[7646]159         END IF
[7698]160         !
161         zarho = SUM( area(:,:) * zbotpres(:,:) )
[7646]162         !                                         
163         IF( lk_mpp )   CALL mpp_sum( zarho )
164         zssh_steric = - zarho / area_tot
165         CALL iom_put( 'sshthster', zssh_steric )
[1756]166     
[7646]167         !                                         ! steric sea surface height
168         CALL eos( tsn, zrhd, zrhop, gdept_n(:,:,:) )                 ! now in situ and potential density
[7698]169!$OMP PARALLEL DO schedule(static) private(jj, ji)
170         DO jj = 1, jpj
171            DO ji = 1, jpi
172               zrhop(ji,jj,jpk) = 0._wp
173            END DO
174         END DO
[7646]175         CALL iom_put( 'rhop', zrhop )
176         !
[7698]177!$OMP PARALLEL
178!$OMP DO schedule(static) private(jj, ji)
179         DO jj = 1, jpj
180            DO ji = 1, jpi
181               zbotpres(ji,jj) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice
182            END DO
183         END DO
[7646]184         DO jk = 1, jpkm1
[7698]185!$OMP DO schedule(static) private(jj, ji)
186            DO jj = 1, jpj
187               DO ji = 1, jpi
188                  zbotpres(ji,jj) = zbotpres(ji,jj) + e3t_n(ji,jj,jk) * zrhd(ji,jj,jk)
189               END DO
190            END DO
[7646]191         END DO
192         IF( ln_linssh ) THEN
193            IF ( ln_isfcav ) THEN
[7698]194!$OMP DO schedule(static) private(jj, ji)
[7646]195               DO ji = 1,jpi
196                  DO jj = 1,jpj
197                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
198                  END DO
[5120]199               END DO
[7646]200            ELSE
[7698]201!$OMP DO schedule(static) private(jj, ji)
202               DO jj = 1, jpj
203                  DO ji = 1, jpi
204                     zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,1)
205                  END DO
206               END DO
[7646]207            END IF
[5120]208         END IF
[7698]209!$OMP END PARALLEL
[7646]210         !   
[7698]211         zarho = SUM( area(:,:) * zbotpres(:,:) )
[7646]212         IF( lk_mpp )   CALL mpp_sum( zarho )
213         zssh_steric = - zarho / area_tot
214         CALL iom_put( 'sshsteric', zssh_steric )
[1756]215     
[7646]216         !                                         ! ocean bottom pressure
217         zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
[7698]218!$OMP PARALLEL DO schedule(static) private(jj, ji)
219         DO jj = 1, jpj
220            DO ji = 1, jpi
221               zbotpres(ji,jj) = zztmp * ( zbotpres(ji,jj) + sshn(ji,jj) + thick0(ji,jj) )
222            END DO
223         END DO
[7646]224         CALL iom_put( 'botpres', zbotpres )
225         !
226      ENDIF
[1756]227
[7646]228      IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' )  .OR. iom_use( 'saltot' )  ) THEN   
229         !                                         ! Mean density anomalie, temperature and salinity
230         ztemp = 0._wp
231         zsal  = 0._wp
232         DO jk = 1, jpkm1
233            DO jj = 1, jpj
234               DO ji = 1, jpi
235                  zztmp = area(ji,jj) * e3t_n(ji,jj,jk)
236                  ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
237                  zsal  = zsal  + zztmp * tsn(ji,jj,jk,jp_sal)
238               END DO
[1756]239            END DO
240         END DO
[7646]241         IF( ln_linssh ) THEN
242            IF( ln_isfcav ) THEN
243               DO ji = 1, jpi
244                  DO jj = 1, jpj
245                     ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) 
246                     zsal  = zsal  + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) 
247                  END DO
[5120]248               END DO
[7646]249            ELSE
250               ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
251               zsal  = zsal  + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
252            END IF
253         ENDIF
254         IF( lk_mpp ) THEN 
255            CALL mpp_sum( ztemp )
256            CALL mpp_sum( zsal  )
[5120]257         END IF
[7646]258         !
259         zmass = rau0 * ( zarho + zvol )                 ! total mass of liquid seawater
260         ztemp = ztemp / zvol                            ! potential temperature in liquid seawater
261         zsal  = zsal  / zvol                            ! Salinity of liquid seawater
262         !
263         CALL iom_put( 'masstot', zmass )
264         CALL iom_put( 'temptot', ztemp )
265         CALL iom_put( 'saltot' , zsal  )
266         !
[1756]267      ENDIF
[7646]268
269      IF( iom_use( 'tnpeo' )) THEN   
270      ! Work done against stratification by vertical mixing
271      ! Exclude points where rn2 is negative as convection kicks in here and
272      ! work is not being done against stratification
273          CALL wrk_alloc( jpi, jpj, zpe )
[7698]274!$OMP PARALLEL DO schedule(static) private(jj,ji)
275          DO jj = 1, jpj
276             DO ji = 1, jpi
277                zpe(ji,jj) = 0._wp
278             END DO
279          END DO
[7646]280          IF( lk_zdfddm ) THEN
[7698]281!$OMP PARALLEL DO schedule(static) private(ji,jj,jk,zrw,zaw,zbw)
[7646]282             DO ji=1,jpi
283                DO jj=1,jpj
284                   DO jk=1,jpk
285                      zrw =   ( gdepw_n(ji,jj,jk  ) - gdept_n(ji,jj,jk) )   &
286                         &  / ( gdept_n(ji,jj,jk-1) - gdept_n(ji,jj,jk) )
287                      !
288                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
289                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
290                      !
291                      zpe(ji, jj) = zpe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * &
292                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  &
293                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
294
295                   ENDDO
296                ENDDO
297             ENDDO
298          ELSE
[7698]299!$OMP PARALLEL DO schedule(static) private(ji,jj,jk)
[7646]300             DO ji = 1, jpi
301                DO jj = 1, jpj
302                   DO jk = 1, jpk
303                       zpe(ji,jj) = zpe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * e3w_n(ji, jj, jk)
304                   ENDDO
305                ENDDO
306             ENDDO
307          ENDIF
308          CALL lbc_lnk( zpe, 'T', 1._wp)         
309          CALL iom_put( 'tnpeo', zpe )
310          CALL wrk_dealloc( jpi, jpj, zpe )
311      ENDIF
[1756]312      !
[7646]313      IF( l_ar5 ) THEN
314        CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres )
315        CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    )
316        CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 )
317      ENDIF
[1756]318      !
[3294]319      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5')
320      !
[1756]321   END SUBROUTINE dia_ar5
322
[7646]323   SUBROUTINE dia_ar5_hst( ktra, cptr, pua, pva ) 
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)   :: pua   ! 3D input array of advection/diffusion
333      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in)   :: pva   ! 3D input array of advection/diffusion
334      !
335      INTEGER    ::  ji, jj, jk
336      REAL(wp), POINTER, DIMENSION(:,:)  :: z2d
[1756]337
[7646]338   
339
340      CALL wrk_alloc( jpi, jpj, z2d )
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( 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( 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       CALL wrk_dealloc( jpi, jpj, z2d )
378
379   END SUBROUTINE dia_ar5_hst
380
381
[1756]382   SUBROUTINE dia_ar5_init
383      !!----------------------------------------------------------------------
384      !!                  ***  ROUTINE dia_ar5_init  ***
385      !!                   
[2528]386      !! ** Purpose :   initialization for AR5 diagnostic computation
[1756]387      !!----------------------------------------------------------------------
388      INTEGER  ::   inum
389      INTEGER  ::   ik
390      INTEGER  ::   ji, jj, jk  ! dummy loop indices
[7698]391      REAL(wp) ::   zztmp, zsum 
[2715]392      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zsaldta   ! Jan/Dec levitus salinity
[5253]393      !
[1756]394      !!----------------------------------------------------------------------
395      !
[3294]396      IF( nn_timing == 1 )   CALL timing_start('dia_ar5_init')
397      !
[7646]398      l_ar5 = .FALSE.
399      IF(   iom_use( 'voltot'  ) .OR. iom_use( 'sshtot'    )  .OR. iom_use( 'sshdyn' )  .OR.  & 
400         &  iom_use( 'masstot' ) .OR. iom_use( 'temptot'   )  .OR. iom_use( 'saltot' ) .OR.  &   
401         &  iom_use( 'botpres' ) .OR. iom_use( 'sshthster' )  .OR. iom_use( 'sshsteric' )  ) L_ar5 = .TRUE.
402 
403      IF( l_ar5 ) THEN
404         !
405         CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )
406         !                                      ! allocate dia_ar5 arrays
407         IF( dia_ar5_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
[2715]408
[7698]409!$OMP PARALLEL DO schedule(static) private(jj, ji)
410         DO jj = 1, jpj
411            DO ji = 1, jpi
412               area(ji,jj) = e1e2t(ji,jj) * tmask_i(ji,jj)
413            END DO
414         END DO
[1756]415
[7646]416         area_tot = SUM( area(:,:) )   ;   IF( lk_mpp )   CALL mpp_sum( area_tot )
[1756]417
[7646]418         vol0        = 0._wp
[7698]419!$OMP PARALLEL
420!$OMP DO schedule(static) private(jj, ji)
421         DO jj = 1, jpj
422            DO ji = 1, jpi
423               thick0(ji,jj) = 0._wp
424            END DO
425         END DO
[7646]426         DO jk = 1, jpkm1
[7698]427!$OMP DO schedule(static) private(jj, ji, zsum)
428            DO jj = 1, jpj
429               DO ji = 1, jpi
430                  zsum = area (ji,jj) * tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
431               END DO
432            END DO
433            vol0        = vol0        + zsum
434!$OMP DO schedule(static) private(jj, ji)
435            DO jj = 1, jpj
436               DO ji = 1, jpi
437                  thick0(ji,jj) = thick0(ji,jj) + tmask_i(ji,jj) * tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
438               END DO
439            END DO
[7646]440         END DO
[7698]441!$OMP END PARALLEL
[7646]442         IF( lk_mpp )   CALL mpp_sum( vol0 )
[5253]443
[7646]444         CALL iom_open ( 'sali_ref_clim_monthly', inum )
445         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1  )
446         CALL iom_get  ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
447         CALL iom_close( inum )
[6665]448
[7698]449!$OMP PARALLEL
450!$OMP DO schedule(static) private(jk, jj, ji)
451         DO jk = 1, jpk
452            DO jj = 1, jpj
453               DO ji = 1, jpi
454                  sn0(ji,jj,jk) = 0.5_wp * ( zsaldta(ji,jj,jk,1) + zsaldta(ji,jj,jk,2) )       
455                  sn0(ji,jj,jk) = sn0(ji,jj,jk) * tmask(ji,jj,jk)
456               END DO
457            END DO
458         END DO
[7646]459         IF( ln_zps ) THEN               ! z-coord. partial steps
[7698]460!$OMP DO schedule(static) private(jj, ji, ik, zztmp)
[7646]461            DO jj = 1, jpj               ! interpolation of salinity at the last ocean level (i.e. the partial step)
462               DO ji = 1, jpi
463                  ik = mbkt(ji,jj)
464                  IF( ik > 1 ) THEN
465                     zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
466                     sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
467                  ENDIF
468               END DO
[1756]469            END DO
[7646]470         ENDIF
[7698]471!$OMP END PARALLEL
[7646]472         !
473         CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )
474         !
[1756]475      ENDIF
476      !
[3294]477      IF( nn_timing == 1 )   CALL timing_stop('dia_ar5_init')
478      !
[1756]479   END SUBROUTINE dia_ar5_init
480
481   !!======================================================================
482END MODULE diaar5
Note: See TracBrowser for help on using the repository browser.