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.
diaptr.F90 in branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/DIA/diaptr.F90 @ 2364

Last change on this file since 2364 was 2364, checked in by acc, 13 years ago

Added basic NetCDF4 chunking and compression support (key_netcdf4). See ticket #754

  • Property svn:keywords set to Id
File size: 47.9 KB
Line 
1MODULE diaptr
2   !!======================================================================
3   !!                       ***  MODULE  diaptr  ***
4   !! Ocean physics:  Computes meridonal transports and zonal means
5   !!=====================================================================
6   !! History :  1.0  ! 2003-09  (C. Talandier, G. Madec)  Original code
7   !!            2.0  ! 2006-01  (A. Biastoch)  Allow sub-basins computation
8   !!            3.2  ! 2003-03  (O. Marti, S. Flavoni) Add fields
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   dia_ptr      : Poleward Transport Diagnostics module
13   !!   dia_ptr_init : Initialization, namelist read
14   !!   dia_ptr_wri  : Output of poleward fluxes
15   !!   ptr_vjk      : "zonal" sum computation of a "meridional" flux array
16   !!   ptr_tjk      : "zonal" mean computation of a tracer field
17   !!   ptr_vj       : "zonal" and vertical sum computation of a "meridional"
18   !!                : flux array; Generic interface: ptr_vj_3d, ptr_vj_2d
19   !!----------------------------------------------------------------------
20   USE oce           ! ocean dynamics and active tracers
21   USE dom_oce       ! ocean space and time domain
22   USE phycst        ! physical constants
23   USE ldftra_oce    ! ocean active tracers: lateral physics
24   USE dianam
25   USE iom
26   USE ioipsl         
27   USE in_out_manager
28   USE lib_mpp
29   USE lbclnk
30
31   IMPLICIT NONE
32   PRIVATE
33
34   INTERFACE ptr_vj
35      MODULE PROCEDURE ptr_vj_3d, ptr_vj_2d
36   END INTERFACE
37
38   PUBLIC   dia_ptr_init   ! call in opa module
39   PUBLIC   dia_ptr        ! call in step module
40   PUBLIC   ptr_vj         ! call by tra_ldf & tra_adv routines
41   PUBLIC   ptr_vjk        ! call by tra_ldf & tra_adv routines
42
43   !                                           !!** namelist  namptr  **
44   LOGICAL , PUBLIC ::   ln_diaptr  = .FALSE.   !: Poleward transport flag (T) or not (F)
45   LOGICAL , PUBLIC ::   ln_subbas  = .FALSE.   !: Atlantic/Pacific/Indian basins calculation
46   LOGICAL , PUBLIC ::   ln_diaznl  = .FALSE.   !: Add zonal means and meridional stream functions
47   LOGICAL , PUBLIC ::   ln_ptrcomp = .FALSE.   !: Add decomposition : overturning (and gyre, soon ...)
48   INTEGER , PUBLIC ::   nf_ptr     = 15        !: frequency of ptr computation
49   INTEGER , PUBLIC ::   nf_ptr_wri = 15        !: frequency of ptr outputs
50
51   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   abasin, pbasin, ibasin, dbasin, sbasin   !: Sub basin masks
52
53   !                                                               !!! poleward heat and salt transport
54   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_adv    , pst_adv       !: advection
55   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ldf    , pst_ldf       !: lateral diffusion
56   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ove_glo, pst_ove_glo   !: global       overturning
57   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ove_atl, pst_ove_atl   !: Atlantic     overturning
58   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ove_pac, pst_ove_pac   !: Pacific      overturning
59   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ove_ind, pst_ove_ind   !: Indian       overturning
60   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_ove_ipc, pst_ove_ipc   !: Indo-Pacific overturning
61   REAL(wp), PUBLIC, DIMENSION(jpj) ::   ht_glo, ht_atl, ht_ind, ht_pac, ht_ipc   !: heat
62   REAL(wp), PUBLIC, DIMENSION(jpj) ::   st_glo, st_atl, st_ind, st_pac, st_ipc   !: salt
63
64   INTEGER ::   niter
65   INTEGER ::   nidom_ptr
66   INTEGER ::   numptr                                              !: logical unit for Poleward TRansports
67
68   REAL(wp), DIMENSION(jpj,jpk) ::   tn_jk_glo  , sn_jk_glo       ! global       i-mean temperature and salinity
69   REAL(wp), DIMENSION(jpj,jpk) ::   tn_jk_atl  , sn_jk_atl       ! Atlantic               -              -
70   REAL(wp), DIMENSION(jpj,jpk) ::   tn_jk_pac  , sn_jk_pac       ! Pacific                -              -
71   REAL(wp), DIMENSION(jpj,jpk) ::   tn_jk_ind  , sn_jk_ind       ! Indian                 -              -
72   REAL(wp), DIMENSION(jpj,jpk) ::   tn_jk_ipc  , sn_jk_ipc       ! Indo-Pacific           -              -
73   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_glo                    ! global       "meridional" Stream-Function
74   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_atl                    ! Atlantic               -              -
75   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_pac                    ! Pacific                -              -
76   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_ind                    ! Indian                 -              -
77   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_ipc                    ! Indo-Pacific           -              -
78   REAL(wp), DIMENSION(jpj,jpk) ::   surf_jk_glo, surf_jk_r_glo   ! surface of global       i-section and its inverse
79   REAL(wp), DIMENSION(jpj,jpk) ::   surf_jk_atl, surf_jk_r_atl   ! surface of Atlantic          -              -
80   REAL(wp), DIMENSION(jpj,jpk) ::   surf_jk_pac, surf_jk_r_pac   ! surface of Pacific           -              -
81   REAL(wp), DIMENSION(jpj,jpk) ::   surf_jk_ind, surf_jk_r_ind   ! surface of Indian            -              -
82   REAL(wp), DIMENSION(jpj,jpk) ::   surf_jk_ipc, surf_jk_r_ipc   ! surface of Indo-Pacific      -              -
83#if defined key_diaeiv
84   !                                                               !!! eddy induced velocity (bolus)
85   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_eiv_glo, pst_eiv_glo   !: global       poleward heat and salt bolus advection
86   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_eiv_atl, pst_eiv_atl   !: Atlantic         -                           -
87   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_eiv_pac, pst_eiv_pac   !: Pacific          -                           -
88   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_eiv_ind, pst_eiv_ind   !: Indian           -                           -
89   REAL(wp), PUBLIC, DIMENSION(jpj) ::   pht_eiv_ipc, pst_eiv_ipc   !: Indo-Pacific     -                           -
90
91   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_eiv_glo   ! global       "meridional" bolus Stream-Function
92   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_eiv_atl   ! Atlantic          -                   -
93   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_eiv_pac   ! Pacific           -                   -
94   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_eiv_ind   ! Indian            -                   -
95   REAL(wp), DIMENSION(jpj,jpk) ::   v_msf_eiv_ipc   ! Indo-Pacific      -                   -
96#endif
97 
98   !! * Substitutions
99#  include "domzgr_substitute.h90"
100#  include "vectopt_loop_substitute.h90"
101   !!----------------------------------------------------------------------
102   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
103   !! $Id$
104   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
105   !!----------------------------------------------------------------------
106
107CONTAINS
108
109   FUNCTION ptr_vj_3d( pva )   RESULT ( p_fval )
110      !!----------------------------------------------------------------------
111      !!                    ***  ROUTINE ptr_vj_3d  ***
112      !!
113      !! ** Purpose :   "zonal" and vertical sum computation of a "meridional"
114      !!              flux array
115      !!
116      !! ** Method  : - i-k sum of pva using the interior 2D vmask (vmask_i).
117      !!              pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v)
118      !!
119      !! ** Action  : - p_fval: i-k-mean poleward flux of pva
120      !!----------------------------------------------------------------------
121      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk) ::   pva   ! mask flux array at V-point
122      !!
123      INTEGER                  ::   ji, jj, jk   ! dummy loop arguments
124      INTEGER                  ::   ijpj         ! ???
125      REAL(wp), DIMENSION(jpj) ::   p_fval       ! function value
126      !!--------------------------------------------------------------------
127      !
128      ijpj = jpj
129      p_fval(:) = 0.e0
130      DO jk = 1, jpkm1
131         DO jj = 2, jpjm1
132            DO ji = fs_2, fs_jpim1   ! Vector opt.
133               p_fval(jj) = p_fval(jj) + pva(ji,jj,jk) * tmask_i(ji,jj) 
134            END DO
135         END DO
136      END DO
137      !
138#if defined key_mpp_mpi
139      CALL mpp_sum( p_fval, ijpj, ncomm_znl)     !!bug  I presume
140#endif
141      !
142   END FUNCTION ptr_vj_3d
143
144
145   FUNCTION ptr_vj_2d( pva )   RESULT ( p_fval )
146      !!----------------------------------------------------------------------
147      !!                    ***  ROUTINE ptr_vj_2d  ***
148      !!
149      !! ** Purpose :   "zonal" and vertical sum computation of a "meridional"
150      !!      flux array
151      !!
152      !! ** Method  : - i-k sum of pva using the interior 2D vmask (vmask_i).
153      !!      pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v)
154      !!
155      !! ** Action  : - p_fval: i-k-mean poleward flux of pva
156      !!----------------------------------------------------------------------
157      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj) ::   pva   ! mask flux array at V-point
158      !!
159      INTEGER                  ::   ji,jj    ! dummy loop arguments
160      INTEGER                  ::   ijpj     ! ???
161      REAL(wp), DIMENSION(jpj) ::   p_fval   ! function value
162      !!--------------------------------------------------------------------
163      !
164      ijpj = jpj
165      p_fval(:) = 0.e0
166      DO jj = 2, jpjm1
167         DO ji = nldi, nlei   ! No vector optimisation here. Better use a mask ?
168            p_fval(jj) = p_fval(jj) + pva(ji,jj) * tmask_i(ji,jj)
169         END DO
170      END DO
171      !
172#if defined key_mpp_mpi
173      CALL mpp_sum( p_fval, ijpj, ncomm_znl )     !!bug  I presume
174#endif
175      !
176   END FUNCTION ptr_vj_2d
177
178
179   FUNCTION ptr_vjk( pva, bmask )   RESULT ( p_fval )
180      !!----------------------------------------------------------------------
181      !!                    ***  ROUTINE ptr_vjk  ***
182      !!
183      !! ** Purpose :   "zonal" sum computation of a "meridional" flux array
184      !!
185      !! ** Method  : - i-sum of pva using the interior 2D vmask (vmask_i).
186      !!      pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v)
187      !!
188      !! ** Action  : - p_fval: i-k-mean poleward flux of pva
189      !!----------------------------------------------------------------------
190      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk)           ::   pva     ! mask flux array at V-point
191      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj)    , OPTIONAL ::   bmask   ! Optional 2D basin mask
192      !!
193      INTEGER                      ::   ji, jj, jk   ! dummy loop arguments
194      REAL(wp), DIMENSION(jpj,jpk) ::   p_fval       ! return function value
195#if defined key_mpp_mpi
196      INTEGER, DIMENSION(1) ::   ish
197      INTEGER, DIMENSION(2) ::   ish2
198      REAL(wp), DIMENSION(jpj*jpk) ::   zwork   ! 1D workspace
199#endif
200      !!--------------------------------------------------------------------
201      !
202      p_fval(:,:) = 0.e0
203      !
204      IF( PRESENT( bmask ) ) THEN
205         DO jk = 1, jpkm1
206            DO jj = 2, jpjm1
207!!gm here, use of tmask_i  ==> no need of loop over nldi, nlei....
208               DO ji =  nldi, nlei   ! No vector optimisation here. Better use a mask ?
209                  p_fval(jj,jk) = p_fval(jj,jk) + pva(ji,jj,jk) * e1v(ji,jj) * fse3v(ji,jj,jk)   &
210                     &                                          * tmask_i(ji,jj) * bmask(ji,jj)
211               END DO
212            END DO
213         END DO
214      ELSE
215         DO jk = 1, jpkm1
216            DO jj = 2, jpjm1
217               DO ji =  nldi, nlei   ! No vector optimisation here. Better use a mask ?
218                  p_fval(jj,jk) = p_fval(jj,jk) + pva(ji,jj,jk) * e1v(ji,jj) * fse3v(ji,jj,jk)   &
219                     &                                          * tmask_i(ji,jj)
220               END DO
221            END DO
222         END DO
223      END IF
224      !
225#if defined key_mpp_mpi
226      ish(1) = jpj*jpk   ;   ish2(1) = jpj   ;   ish2(2) = jpk
227      zwork(:) = RESHAPE( p_fval, ish )
228      CALL mpp_sum( zwork, jpj*jpk, ncomm_znl )
229      p_fval(:,:) = RESHAPE( zwork, ish2 )
230#endif
231      !
232   END FUNCTION ptr_vjk
233
234
235   FUNCTION ptr_tjk( pta, bmask )   RESULT ( p_fval )
236      !!----------------------------------------------------------------------
237      !!                    ***  ROUTINE ptr_tjk  ***
238      !!
239      !! ** Purpose :   "zonal" mean computation of a tracer field
240      !!
241      !! ** Method  : - i-sum of mj(pta) using tmask
242      !!      multiplied by the inverse of the surface of the "zonal" ocean
243      !!      section
244      !!
245      !! ** Action  : - p_fval: i-k-mean poleward flux of pta
246      !!----------------------------------------------------------------------
247      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk) ::   pta    ! tracer flux array at T-point
248      REAL(wp) , INTENT(in), DIMENSION(jpi,jpj), OPTIONAL :: bmask ! Optional 2D basin mask
249      !!
250      INTEGER                     ::   ji, jj, jk   ! dummy loop arguments
251      REAL(wp),DIMENSION(jpj,jpk) ::   p_fval       ! return function value
252#if defined key_mpp_mpi
253      INTEGER, DIMENSION(1) ::   ish
254      INTEGER, DIMENSION(2) ::   ish2
255      REAL(wp),DIMENSION(jpj*jpk) ::   zwork   ! 1D workspace
256#endif
257      !!--------------------------------------------------------------------
258      !
259      p_fval(:,:) = 0.e0
260      IF (PRESENT (bmask)) THEN
261         DO jk = 1, jpkm1
262            DO jj = 2, jpjm1
263               DO ji =  nldi, nlei   ! No vector optimisation here. Better use a mask ?
264                  p_fval(jj,jk) = p_fval(jj,jk) + pta(ji,jj,jk)                  &
265                     &                          * e1t(ji,jj) * fse3t(ji,jj,jk)   &
266                     &                          * tmask_i(ji,jj)                 &
267                     &                          * bmask(ji,jj)
268               END DO
269            END DO
270         END DO
271      ELSE
272         DO jk = 1, jpkm1
273            DO jj = 2, jpjm1
274               DO ji =  nldi, nlei   ! No vector optimisation here. Better use a mask ?
275                  p_fval(jj,jk) = p_fval(jj,jk) + pta(ji,jj,jk)                  &
276                     &                          * e1t(ji,jj) * fse3t(ji,jj,jk)   &
277                     &                          * tmask_i(ji,jj)
278               END DO
279            END DO
280         END DO
281      END IF
282      p_fval(:,:) = p_fval(:,:) * 0.5
283#if defined key_mpp_mpi
284      ish(1) = jpj*jpk   ;   ish2(1) = jpj   ;   ish2(2) = jpk
285      zwork(:)= RESHAPE( p_fval, ish )
286      CALL mpp_sum( zwork, jpj*jpk, ncomm_znl )
287      p_fval(:,:)= RESHAPE(zwork,ish2)
288#endif
289      !
290   END FUNCTION ptr_tjk
291
292
293   SUBROUTINE dia_ptr( kt )
294      !!----------------------------------------------------------------------
295      !!                  ***  ROUTINE dia_ptr  ***
296      !!----------------------------------------------------------------------
297      INTEGER, INTENT(in) ::   kt   ! ocean time step index
298      !!
299      INTEGER  ::   jk, jj, ji   ! dummy loop
300      REAL(wp) ::   zsverdrup    ! conversion from m3/s to Sverdrup
301      REAL(wp) ::   zpwatt       ! conversion from W    to PW
302      REAL(wp) ::   zggram       ! conversion from g    to Pg
303      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   vt, vs   ! 3D workspace
304      !!----------------------------------------------------------------------
305
306      IF( kt == nit000 .OR. MOD( kt, nf_ptr ) == 0 )   THEN
307
308         IF ( MOD( kt, nf_ptr ) == 0 ) THEN
309
310            zsverdrup = 1.e-6
311            zpwatt    = 1.e-15
312            zggram    = 1.e-6
313
314            IF ( ln_diaznl ) THEN
315               ! "zonal" mean temperature and salinity at V-points
316               tn_jk_glo(:,:) = ptr_tjk( tn(:,:,:) ) * surf_jk_r_glo(:,:)
317               sn_jk_glo(:,:) = ptr_tjk( sn(:,:,:) ) * surf_jk_r_glo(:,:)
318
319               IF (ln_subbas) THEN
320                  tn_jk_atl(:,:) = ptr_tjk( tn(:,:,:), abasin(:,:) ) * surf_jk_r_atl(:,:)
321                  sn_jk_atl(:,:) = ptr_tjk( sn(:,:,:), abasin(:,:) ) * surf_jk_r_atl(:,:)
322                  tn_jk_pac(:,:) = ptr_tjk( tn(:,:,:), pbasin(:,:) ) * surf_jk_r_pac(:,:)
323                  sn_jk_pac(:,:) = ptr_tjk( sn(:,:,:), pbasin(:,:) ) * surf_jk_r_pac(:,:)
324                  tn_jk_ind(:,:) = ptr_tjk( tn(:,:,:), ibasin(:,:) ) * surf_jk_r_ind(:,:)
325                  sn_jk_ind(:,:) = ptr_tjk( sn(:,:,:), ibasin(:,:) ) * surf_jk_r_ind(:,:)
326                  tn_jk_ipc(:,:) = ptr_tjk( tn(:,:,:), dbasin(:,:) ) * surf_jk_r_ipc(:,:)
327                  sn_jk_ipc(:,:) = ptr_tjk( sn(:,:,:), dbasin(:,:) ) * surf_jk_r_ipc(:,:)
328               ENDIF
329            ENDIF
330
331            !--------------------------------------------------------
332            ! overturning calculation:
333
334            ! horizontal integral and vertical dz
335
336#if defined key_diaeiv
337            v_msf_glo(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:) ) 
338            IF( ln_subbas .AND. ln_diaznl ) THEN
339               v_msf_atl(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:), abasin(:,:)*sbasin(:,:) ) 
340               v_msf_pac(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:), pbasin(:,:)*sbasin(:,:) ) 
341               v_msf_ind(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:), ibasin(:,:)*sbasin(:,:) ) 
342               v_msf_ipc(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:), dbasin(:,:)*sbasin(:,:) ) 
343            ENDIF
344#else
345            v_msf_glo(:,:) = ptr_vjk( vn(:,:,:) ) 
346            IF( ln_subbas .AND. ln_diaznl ) THEN
347               v_msf_atl(:,:) = ptr_vjk( vn(:,:,:), abasin(:,:)*sbasin(:,:) ) 
348               v_msf_pac(:,:) = ptr_vjk( vn(:,:,:), pbasin(:,:)*sbasin(:,:) ) 
349               v_msf_ind(:,:) = ptr_vjk( vn(:,:,:), ibasin(:,:)*sbasin(:,:) ) 
350               v_msf_ipc(:,:) = ptr_vjk( vn(:,:,:), dbasin(:,:)*sbasin(:,:) ) 
351            ENDIF
352#endif
353
354#if defined key_diaeiv
355            v_msf_eiv_glo(:,:) = ptr_vjk( v_eiv(:,:,:) )
356            IF (ln_subbas ) THEN
357               v_msf_eiv_atl(:,:) = ptr_vjk( v_eiv(:,:,:), abasin(:,:)*sbasin(:,:) )
358               v_msf_eiv_pac(:,:) = ptr_vjk( v_eiv(:,:,:), pbasin(:,:)*sbasin(:,:) )
359               v_msf_eiv_ind(:,:) = ptr_vjk( v_eiv(:,:,:), ibasin(:,:)*sbasin(:,:) )
360               v_msf_eiv_ipc(:,:) = ptr_vjk( v_eiv(:,:,:), dbasin(:,:)*sbasin(:,:) )
361            END IF
362#endif
363
364            ! Transports
365            ! T times V on T points (include bolus velocities)
366#if defined key_diaeiv 
367            DO jj = 2, jpj
368               DO ji = 1, jpi
369                  vt(ji,jj,:) = tn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) + u_eiv(ji,jj,:) + u_eiv(ji,jj-1,:) )*0.5
370                  vs(ji,jj,:) = sn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) + v_eiv(ji,jj,:) + v_eiv(ji,jj-1,:) )*0.5
371               END DO
372            END DO
373#else
374            DO jj = 2, jpj
375               DO ji = 1, jpi
376                  vt(ji,jj,:) = tn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) )*0.5
377                  vs(ji,jj,:) = sn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) )*0.5
378               END DO
379            END DO
380#endif
381            CALL lbc_lnk( vs, 'V', -1. )   ;   CALL lbc_lnk( vt, 'V', -1. )
382
383            ht_glo(:) = SUM( ptr_vjk( vt(:,:,:)), 2 )
384            st_glo(:) = SUM( ptr_vjk( vs(:,:,:)), 2 )
385
386            IF ( ln_subbas ) THEN
387               ht_atl(:) = SUM( ptr_vjk( vt (:,:,:), abasin(:,:)*sbasin(:,:)), 2 )
388               ht_pac(:) = SUM( ptr_vjk( vt (:,:,:), pbasin(:,:)*sbasin(:,:)), 2 )
389               ht_ind(:) = SUM( ptr_vjk( vt (:,:,:), ibasin(:,:)*sbasin(:,:)), 2 )
390               ht_ipc(:) = SUM( ptr_vjk( vt (:,:,:), dbasin(:,:)*sbasin(:,:)), 2 )
391               st_atl(:) = SUM( ptr_vjk( vs (:,:,:), abasin(:,:)*sbasin(:,:)), 2 )
392               st_pac(:) = SUM( ptr_vjk( vs (:,:,:), pbasin(:,:)*sbasin(:,:)), 2 )
393               st_ind(:) = SUM( ptr_vjk( vs (:,:,:), ibasin(:,:)*sbasin(:,:)), 2 )
394               st_ipc(:) = SUM( ptr_vjk( vs (:,:,:), dbasin(:,:)*sbasin(:,:)), 2 )
395            ENDIF
396
397            ! poleward tracer transports:
398            ! overturning components:
399            IF ( ln_ptrcomp ) THEN
400               pht_ove_glo(:) = SUM( v_msf_glo(:,:) * tn_jk_glo(:,:), 2 )   ! SUM over jk
401               pst_ove_glo(:) = SUM( v_msf_glo(:,:) * sn_jk_glo(:,:), 2 ) 
402               IF ( ln_subbas ) THEN
403                  pht_ove_atl(:) = SUM( v_msf_atl(:,:) * tn_jk_atl(:,:), 2 )   ! SUM over jk
404                  pst_ove_atl(:) = SUM( v_msf_atl(:,:) * sn_jk_atl(:,:), 2 ) 
405                  pht_ove_pac(:) = SUM( v_msf_pac(:,:) * tn_jk_pac(:,:), 2 )   ! SUM over jk
406                  pst_ove_pac(:) = SUM( v_msf_pac(:,:) * sn_jk_pac(:,:), 2 ) 
407                  pht_ove_ind(:) = SUM( v_msf_ind(:,:) * tn_jk_ind(:,:), 2 )   ! SUM over jk
408                  pst_ove_ind(:) = SUM( v_msf_ind(:,:) * sn_jk_ind(:,:), 2 ) 
409                  pht_ove_ipc(:) = SUM( v_msf_ipc(:,:) * tn_jk_ipc(:,:), 2 )   ! SUM over jk
410                  pst_ove_ipc(:) = SUM( v_msf_ipc(:,:) * sn_jk_ipc(:,:), 2 ) 
411               END IF
412            END IF
413
414            ! Bolus component
415#if defined key_diaeiv
416            pht_eiv_glo(:) = SUM( v_msf_eiv_glo(:,:) * tn_jk_glo(:,:), 2 )   ! SUM over jk
417            pst_eiv_glo(:) = SUM( v_msf_eiv_glo(:,:) * sn_jk_glo(:,:), 2 )   ! SUM over jk
418            IF ( ln_subbas ) THEN
419               pht_eiv_atl(:) = SUM( v_msf_eiv_glo(:,:) * tn_jk_atl(:,:), 2 )   ! SUM over jk
420               pst_eiv_atl(:) = SUM( v_msf_eiv_glo(:,:) * sn_jk_atl(:,:), 2 )   ! SUM over jk
421               pht_eiv_pac(:) = SUM( v_msf_eiv_pac(:,:) * tn_jk_pac(:,:), 2 )   ! SUM over jk
422               pst_eiv_pac(:) = SUM( v_msf_eiv_pac(:,:) * sn_jk_pac(:,:), 2 )   ! SUM over jk
423               pht_eiv_ind(:) = SUM( v_msf_eiv_ind(:,:) * tn_jk_ind(:,:), 2 )   ! SUM over jk
424               pst_eiv_ind(:) = SUM( v_msf_eiv_ind(:,:) * sn_jk_ind(:,:), 2 )   ! SUM over jk
425               pht_eiv_ipc(:) = SUM( v_msf_eiv_ipc(:,:) * tn_jk_ipc(:,:), 2 )   ! SUM over jk
426               pst_eiv_ipc(:) = SUM( v_msf_eiv_ipc(:,:) * sn_jk_ipc(:,:), 2 )   ! SUM over jk
427            ENDIF
428#endif
429
430            ! conversion in PW and G g
431            zpwatt = zpwatt * rau0 * rcp
432            pht_adv(:) = pht_adv(:) * zpwatt 
433            pht_ldf(:) = pht_ldf(:) * zpwatt
434            pst_adv(:) = pst_adv(:) * zggram
435            pst_ldf(:) = pst_ldf(:) * zggram
436            IF ( ln_ptrcomp ) THEN
437               pht_ove_glo(:) = pht_ove_glo(:) * zpwatt
438               pst_ove_glo(:) = pst_ove_glo(:) * zggram
439            END IF
440#if defined key_diaeiv
441            pht_eiv_glo(:) = pht_eiv_glo(:) * zpwatt
442            pst_eiv_glo(:) = pst_eiv_glo(:) * zggram
443#endif
444            IF( ln_subbas ) THEN
445               ht_atl(:) = ht_atl(:) * zpwatt
446               ht_pac(:) = ht_pac(:) * zpwatt
447               ht_ind(:) = ht_ind(:) * zpwatt
448               ht_ipc(:) = ht_ipc(:) * zpwatt
449               st_atl(:) = st_atl(:) * zggram 
450               st_pac(:) = st_pac(:) * zggram
451               st_ind(:) = st_ind(:) * zggram
452               st_ipc(:) = st_ipc(:) * zggram
453            ENDIF
454
455            ! "Meridional" Stream-Function
456            DO jk = 2,jpk 
457               v_msf_glo(:,jk) = v_msf_glo(:,jk-1) + v_msf_glo(:,jk)
458            END DO
459            v_msf_glo(:,:) = v_msf_glo(:,:) * zsverdrup
460#if defined key_diaeiv
461            ! Bolus "Meridional" Stream-Function
462            DO jk = 2,jpk
463               v_msf_eiv_glo(:,jk) = v_msf_eiv_glo(:,jk-1) + v_msf_eiv_glo(:,jk)
464            END DO
465            v_msf_eiv_glo(:,:) = v_msf_eiv_glo(:,:) * zsverdrup
466            IF ( ln_subbas ) THEN
467               DO jk = 2,jpk
468                  v_msf_eiv_atl(:,jk) = v_msf_eiv_atl(:,jk-1) + v_msf_eiv_atl(:,jk)
469                  v_msf_eiv_pac(:,jk) = v_msf_eiv_pac(:,jk-1) + v_msf_eiv_pac(:,jk)
470                  v_msf_eiv_ind(:,jk) = v_msf_eiv_ind(:,jk-1) + v_msf_eiv_ind(:,jk)
471                  v_msf_eiv_ipc(:,jk) = v_msf_eiv_ipc(:,jk-1) + v_msf_eiv_ipc(:,jk)
472               END DO
473            ENDIF
474#endif
475            !
476            IF( ln_subbas .AND. ln_diaznl ) THEN
477               DO jk = 2,jpk 
478                  v_msf_atl(:,jk) = v_msf_atl(:,jk-1) + v_msf_atl(:,jk)
479                  v_msf_pac(:,jk) = v_msf_pac(:,jk-1) + v_msf_pac(:,jk)
480                  v_msf_ind(:,jk) = v_msf_ind(:,jk-1) + v_msf_ind(:,jk)
481                  v_msf_ipc(:,jk) = v_msf_ipc(:,jk-1) + v_msf_ipc(:,jk)
482               END DO
483               v_msf_atl(:,:) = v_msf_atl(:,:) * zsverdrup
484               v_msf_pac(:,:) = v_msf_pac(:,:) * zsverdrup
485               v_msf_ind(:,:) = v_msf_ind(:,:) * zsverdrup
486               v_msf_ipc(:,:) = v_msf_ipc(:,:) * zsverdrup
487            ENDIF
488         ENDIF
489         !
490         CALL dia_ptr_wri( kt )                        ! outputs
491         !
492      ENDIF
493      !
494      IF( kt == nitend )   CALL histclo( numptr )      ! Close the file
495      !
496   END SUBROUTINE dia_ptr
497
498
499   SUBROUTINE dia_ptr_init
500      !!----------------------------------------------------------------------
501      !!                  ***  ROUTINE dia_ptr_init  ***
502      !!                   
503      !! ** Purpose :   Initialization, namelist read
504      !!----------------------------------------------------------------------
505      INTEGER ::   inum       ! temporary logical unit
506#if defined key_mpp_mpi
507      INTEGER, DIMENSION(1) :: iglo, iloc, iabsf, iabsl, ihals, ihale, idid
508#endif
509      !!
510      NAMELIST/namptr/ ln_diaptr, ln_diaznl, ln_subbas, ln_ptrcomp, nf_ptr, nf_ptr_wri
511      !!----------------------------------------------------------------------
512
513      REWIND ( numnam )              ! Read Namelist namptr : poleward transport parameters
514      READ   ( numnam, namptr )
515
516      IF(lwp) THEN                   ! Control print
517         WRITE(numout,*)
518         WRITE(numout,*) 'dia_ptr_init : poleward transport and msf initialization'
519         WRITE(numout,*) '~~~~~~~~~~~~'
520         WRITE(numout,*) '   Namelist namptr : set ptr parameters'
521         WRITE(numout,*) '      Switch for ptr diagnostic (T) or not (F)  ln_diaptr  = ', ln_diaptr
522         WRITE(numout,*) '      Atl/Pac/Ind basins computation            ln_subbas  = ', ln_subbas
523         WRITE(numout,*) '      Frequency of computation                  nf_ptr     = ', nf_ptr
524         WRITE(numout,*) '      Frequency of outputs                      nf_ptr_wri = ', nf_ptr_wri
525      ENDIF
526
527      IF( .NOT. ln_diaptr )   RETURN
528     
529      IF( lk_mpp )   CALL mpp_ini_znl      ! Define MPI communicator for zonal sum
530
531      IF( ln_subbas ) THEN                 ! load sub-basin mask
532         CALL iom_open( 'subbasins', inum )
533         CALL iom_get( inum, jpdom_data, 'atlmsk', abasin )      ! Atlantic basin
534         CALL iom_get( inum, jpdom_data, 'pacmsk', pbasin )      ! Pacific basin
535         CALL iom_get( inum, jpdom_data, 'indmsk', ibasin )      ! Indian basin
536         CALL iom_close( inum )
537         dbasin(:,:) = MAX ( pbasin(:,:), ibasin(:,:) )
538         sbasin(:,:) = tmask (:,:,1)
539         WHERE ( gphit (:,:) < -30.e0) sbasin(:,:) = 0.e0
540      ENDIF
541     
542!!gm CAUTION : this is only valid in fixed volume case !
543
544      ! inverse of the ocean "zonal" v-point section
545      surf_jk_glo(:,:) = ptr_tjk( tmask(:,:,:) )
546      surf_jk_r_glo(:,:) = 0.e0
547      WHERE( surf_jk_glo(:,:) /= 0.e0 )   surf_jk_r_glo(:,:) = 1.e0 / surf_jk_glo(:,:)
548     
549      IF (ln_subbas) THEN
550         surf_jk_atl(:,:) = ptr_tjk( tmask (:,:,:), abasin(:,:) )
551         surf_jk_r_atl(:,:) = 0.e0
552         WHERE( surf_jk_atl(:,:) /= 0.e0 )   surf_jk_r_atl(:,:) = 1.e0 / surf_jk_atl(:,:)
553         !
554         surf_jk_pac(:,:) = ptr_tjk( tmask (:,:,:), pbasin(:,:) )
555         surf_jk_r_pac(:,:) = 0.e0
556         WHERE( surf_jk_pac(:,:) /= 0.e0 )   surf_jk_r_pac(:,:) = 1.e0 / surf_jk_pac(:,:)
557         !
558         surf_jk_ind(:,:) = ptr_tjk( tmask (:,:,:), ibasin(:,:) )
559         surf_jk_r_ind(:,:) = 0.e0
560         WHERE( surf_jk_ind(:,:) /= 0.e0 )   surf_jk_r_ind(:,:) = 1.e0 / surf_jk_ind(:,:)
561         !
562         surf_jk_ipc(:,:) = ptr_tjk( tmask (:,:,:), dbasin(:,:) )
563         surf_jk_r_ipc(:,:) = 0.e0
564         WHERE( surf_jk_ipc(:,:) /= 0.e0 )   surf_jk_r_ipc(:,:) = 1.e0 / surf_jk_ipc(:,:)
565      END IF
566
567     
568      !!----------------------------------------------------------------------
569
570#if defined key_mpp_mpi 
571      iglo (1) = jpjglo
572      iloc (1) = nlcj
573      iabsf(1) = njmppt(narea)
574      iabsl(:) = iabsf(:) + iloc(:) - 1
575      ihals(1) = nldj - 1
576      ihale(1) = nlcj - nlej
577      idid (1) = 2
578
579!-$$      IF(lwp) THEN
580!-$$          WRITE(numout,*)
581!-$$          WRITE(numout,*) 'dia_ptr_init :   iloc  = ', iloc
582!-$$          WRITE(numout,*) '~~~~~~~~~~~~     iabsf = ', iabsf
583!-$$          WRITE(numout,*) '                 ihals = ', ihals
584!-$$          WRITE(numout,*) '                 ihale = ', ihale
585!-$$      ENDIF
586
587      CALL flio_dom_set ( jpnj, nproc/jpni, idid, iglo, iloc, iabsf, iabsl, ihals, ihale, 'BOX', nidom_ptr)
588#else
589      nidom_ptr = FLIO_DOM_NONE
590#endif
591      !
592   END SUBROUTINE dia_ptr_init
593
594
595   SUBROUTINE dia_ptr_wri( kt )
596      !!---------------------------------------------------------------------
597      !!                ***  ROUTINE dia_ptr_wri  ***
598      !!
599      !! ** Purpose :   output of poleward fluxes
600      !!
601      !! ** Method  :   NetCDF file
602      !!----------------------------------------------------------------------
603      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
604      !!
605      INTEGER, SAVE ::   nhoridz, ndepidzt, ndepidzw
606      INTEGER, SAVE :: ndim  , ndim_atl     , ndim_pac     , ndim_ind     , ndim_ipc
607      INTEGER, SAVE ::         ndim_atl_30  , ndim_pac_30  , ndim_ind_30  , ndim_ipc_30
608      INTEGER, SAVE :: ndim_h, ndim_h_atl_30, ndim_h_pac_30, ndim_h_ind_30, ndim_h_ipc_30
609      INTEGER, SAVE, DIMENSION (jpj*jpk) :: ndex  , ndex_atl     , ndex_pac     , ndex_ind     , ndex_ipc
610      INTEGER, SAVE, DIMENSION (jpj*jpk) ::         ndex_atl_30  , ndex_pac_30  , ndex_ind_30  , ndex_ipc_30
611      INTEGER, SAVE, DIMENSION (jpj)     :: ndex_h, ndex_h_atl_30, ndex_h_pac_30, ndex_h_ind_30, ndex_h_ipc_30
612
613      CHARACTER (len=40)       ::   clhstnam, clop, clop_once, cl_comment   ! temporary names
614      INTEGER                  ::   iline, it, itmod, ji, jj, jk            !
615#if defined key_iomput
616      INTEGER                  ::   inum                                    ! temporary logical unit
617#endif
618      REAL(wp)                 ::   zsto, zout, zdt, zjulian                ! temporary scalars
619      REAL(wp), DIMENSION(jpj) ::   zphi, zfoo
620      REAL(wp), DIMENSION(jpj,jpk) :: z_1
621      !!----------------------------------------------------------------------
622
623      ! define time axis
624      it    = kt / nf_ptr
625      itmod = kt - nit000 + 1
626     
627!-$$      IF(lwp) THEN
628!-$$         WRITE(numout,*)
629!-$$         WRITE(numout,*) 'dia_ptr_wri : kt = ', kt, 'it = ', it, ' itmod = ', itmod, ' niter = ', niter
630!-$$         WRITE(numout,*) '~~~~~~~~~~~~'
631!-$$      ENDIF
632
633      ! Initialization
634      ! --------------
635      IF( kt == nit000 ) THEN
636
637         niter = (nit000 - 1) / nf_ptr
638
639!-$$         IF(lwp) THEN
640!-$$            WRITE(numout,*)
641!-$$            WRITE(numout,*) 'dia_ptr_wri : poleward transport and msf writing: initialization , niter = ', niter
642!-$$            WRITE(numout,*) '~~~~~~~~~~~~'
643!-$$         ENDIF
644
645         zdt = rdt
646         IF( nacc == 1 )   zdt = rdtmin
647
648         ! Reference latitude
649         ! ------------------
650         !                                           ! =======================
651         IF( cp_cfg == "orca" ) THEN                 !   ORCA configurations
652            !                                        ! =======================
653
654            IF( jp_cfg == 05  )   iline = 192   ! i-line that passes near the North Pole
655            IF( jp_cfg == 025 )   iline = 384   ! i-line that passes near the North Pole
656            IF( jp_cfg == 1   )   iline =  96   ! i-line that passes near the North Pole
657            IF( jp_cfg == 2   )   iline =  48   ! i-line that passes near the North Pole
658            IF( jp_cfg == 4   )   iline =  24   ! i-line that passes near the North Pole
659            zphi(:) = 0.e0
660            DO ji = mi0(iline), mi1(iline) 
661               zphi(:) = gphiv(ji,:)         ! if iline is in the local domain
662               ! Correct highest latitude for some configurations - will work if domain is parallelized in J ?
663               IF( jp_cfg == 05 ) THEN
664                  DO jj = mj0(jpjdta), mj1(jpjdta) 
665                     zphi( jj ) = zphi(mj0(jpjdta-1)) + (zphi(mj0(jpjdta-1))-zphi(mj0(jpjdta-2)))/2.
666                     zphi( jj ) = MIN( zphi(jj), 90.)
667                  END DO
668               END IF
669               IF( jp_cfg == 1 .OR. jp_cfg == 2 .OR. jp_cfg == 4 ) THEN
670                  DO jj = mj0(jpjdta-1), mj1(jpjdta-1) 
671                     zphi( jj ) = 88.5e0
672                  END DO
673                  DO jj = mj0(jpjdta  ), mj1(jpjdta  ) 
674                     zphi( jj ) = 89.5e0
675                  END DO
676               END IF
677            END DO
678            ! provide the correct zphi to all local domains
679#if defined key_mpp_mpi
680            CALL mpp_sum( zphi, jpj, ncomm_znl )       
681#endif
682
683            !                                        ! =======================
684         ELSE                                        !   OTHER configurations
685            !                                        ! =======================
686            zphi(:) = gphiv(1,:)             ! assume lat/lon coordinate, select the first i-line
687            !
688         ENDIF
689         !
690         ! Work only on westmost processor (will not work if mppini2 is used)
691#if defined key_mpp_mpi
692         IF ( l_znl_root ) THEN 
693#endif
694            !
695            ! OPEN netcdf file
696            ! ----------------
697            ! Define frequency of output and means
698            zsto = nf_ptr * zdt
699            IF( ln_mskland )   THEN    ! put 1.e+20 on land (very expensive!!)
700               clop      = "ave(only(x))"
701               clop_once = "once(only(x))"
702            ELSE                       ! no use of the mask value (require less cpu time)
703               clop      = "ave(x)"       
704               clop_once = "once"
705            ENDIF
706
707            zout = nf_ptr_wri * zdt
708            zfoo(:) = 0.e0
709
710            ! Compute julian date from starting date of the run
711
712            CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian )
713            zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment
714
715#if defined key_iomput
716            ! Requested by IPSL people, use by their postpro...
717            IF(lwp) THEN
718               CALL dia_nam( clhstnam, nf_ptr_wri,' ' )
719               CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, narea )
720               WRITE(inum,*) clhstnam
721               CLOSE(inum)
722            ENDIF
723#endif
724
725            CALL dia_nam( clhstnam, nf_ptr_wri, 'diaptr' )
726            IF(lwp)WRITE( numout,*)" Name of diaptr NETCDF file : ", clhstnam
727
728            ! Horizontal grid : zphi()
729            CALL histbeg(clhstnam, 1, zfoo, jpj, zphi,   &
730               1, 1, 1, jpj, niter, zjulian, zdt*nf_ptr, nhoridz, numptr, domain_id=nidom_ptr, snc4chunks=snc4set)
731            ! Vertical grids : gdept_0, gdepw_0
732            CALL histvert( numptr, "deptht", "Vertical T levels",   &
733               "m", jpk, gdept_0, ndepidzt, "down" )
734            CALL histvert( numptr, "depthw", "Vertical W levels",   &
735               "m", jpk, gdepw_0, ndepidzw, "down" )
736
737            !
738            CALL wheneq ( jpj*jpk, MIN(surf_jk_glo(:,:), 1.e0), 1, 1., ndex  , ndim  )      ! Lat-Depth
739            CALL wheneq ( jpj    , MIN(surf_jk_glo(:,1), 1.e0), 1, 1., ndex_h, ndim_h )     ! Lat
740
741            IF (ln_subbas) THEN
742               z_1 (:,1) = 1.0e0
743               WHERE ( gphit (jpi/2,:) .LT. -30 ) z_1 (:,1) = 0.e0
744               DO jk = 2, jpk
745                  z_1 (:,jk) = z_1 (:,1)
746               END DO
747
748               CALL wheneq ( jpj*jpk, MIN(surf_jk_atl(:,:)         , 1.e0), 1, 1., ndex_atl     , ndim_atl      ) ! Lat-Depth
749               CALL wheneq ( jpj*jpk, MIN(surf_jk_atl(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_atl_30  , ndim_atl_30   ) ! Lat-Depth
750               CALL wheneq ( jpj    , MIN(surf_jk_atl(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_atl_30, ndim_h_atl_30 ) ! Lat
751
752               CALL wheneq ( jpj*jpk, MIN(surf_jk_pac(:,:)         , 1.e0), 1, 1., ndex_pac     , ndim_pac      ) ! Lat-Depth
753               CALL wheneq ( jpj*jpk, MIN(surf_jk_pac(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_pac_30  , ndim_pac_30   ) ! Lat-Depth
754               CALL wheneq ( jpj    , MIN(surf_jk_pac(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_pac_30, ndim_h_pac_30 ) ! Lat
755
756               CALL wheneq ( jpj*jpk, MIN(surf_jk_ind(:,:)         , 1.e0), 1, 1., ndex_ind     , ndim_ind      ) ! Lat-Depth
757               CALL wheneq ( jpj*jpk, MIN(surf_jk_ind(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_ind_30  , ndim_ind_30   ) ! Lat-Depth
758               CALL wheneq ( jpj    , MIN(surf_jk_ind(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_ind_30, ndim_h_ind_30 ) ! Lat
759
760               CALL wheneq ( jpj*jpk, MIN(surf_jk_ipc(:,:)         , 1.e0), 1, 1., ndex_ipc     , ndim_ipc      ) ! Lat-Depth
761               CALL wheneq ( jpj*jpk, MIN(surf_jk_ipc(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_ipc_30  , ndim_ipc_30   ) ! Lat-Depth
762               CALL wheneq ( jpj    , MIN(surf_jk_ipc(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_ipc_30, ndim_h_ipc_30 ) ! Lat
763
764            ENDIF
765
766            !
767#if defined key_diaeiv
768            cl_comment = ' (Bolus part included)'
769#else
770            cl_comment = '                      '
771#endif
772            !  Zonal mean T and S
773
774            IF ( ln_diaznl ) THEN
775               CALL histdef( numptr, "zotemglo", "Zonal Mean Temperature","C" ,   &
776                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
777               CALL histdef( numptr, "zosalglo", "Zonal Mean Salinity","PSU"  ,   &
778                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
779
780               CALL histdef( numptr, "zosrfglo", "Zonal Mean Surface","m^2"   ,   &
781                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout )
782
783               IF (ln_subbas) THEN
784                  CALL histdef( numptr, "zotematl", "Zonal Mean Temperature: Atlantic","C" ,   &
785                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
786                  CALL histdef( numptr, "zosalatl", "Zonal Mean Salinity: Atlantic","PSU"  ,   &
787                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
788                  CALL histdef( numptr, "zosrfatl", "Zonal Mean Surface: Atlantic","m^2"   ,   &
789                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout )
790
791                  CALL histdef( numptr, "zotempac", "Zonal Mean Temperature: Pacific","C"  ,   &
792                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
793                  CALL histdef( numptr, "zosalpac", "Zonal Mean Salinity: Pacific","PSU"   ,   &
794                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
795                  CALL histdef( numptr, "zosrfpac", "Zonal Mean Surface: Pacific","m^2"    ,   &
796                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout )
797
798                  CALL histdef( numptr, "zotemind", "Zonal Mean Temperature: Indian","C"   ,   &
799                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
800                  CALL histdef( numptr, "zosalind", "Zonal Mean Salinity: Indian","PSU"    ,   &
801                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
802                  CALL histdef( numptr, "zosrfind", "Zonal Mean Surface: Indian","m^2"     ,   &
803                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout )
804
805                  CALL histdef( numptr, "zotemipc", "Zonal Mean Temperature: Pacific+Indian","C" ,   &
806                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
807                  CALL histdef( numptr, "zosalipc", "Zonal Mean Salinity: Pacific+Indian","PSU"  ,   &
808                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout )
809                  CALL histdef( numptr, "zosrfipc", "Zonal Mean Surface: Pacific+Indian","m^2"   ,   &
810                     1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout )
811               ENDIF
812
813            ENDIF
814
815            !  Meridional Stream-Function (Eulerian and Bolus)
816
817            CALL histdef( numptr, "zomsfglo", "Meridional Stream-Function: Global"//TRIM(cl_comment),"Sv" ,   &
818               1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
819            IF( ln_subbas .AND. ln_diaznl ) THEN
820               CALL histdef( numptr, "zomsfatl", "Meridional Stream-Function: Atlantic"//TRIM(cl_comment),"Sv" ,   &
821                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
822               CALL histdef( numptr, "zomsfpac", "Meridional Stream-Function: Pacific"//TRIM(cl_comment),"Sv"  ,   &
823                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
824               CALL histdef( numptr, "zomsfind", "Meridional Stream-Function: Indian"//TRIM(cl_comment),"Sv"   ,   &
825                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
826               CALL histdef( numptr, "zomsfipc", "Meridional Stream-Function: Indo-Pacific"//TRIM(cl_comment),"Sv" ,&
827                  1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
828            ENDIF
829
830            !  Heat transport
831
832            CALL histdef( numptr, "sophtadv", "Advective Heat Transport"      ,   &
833               "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
834            CALL histdef( numptr, "sophtldf", "Diffusive Heat Transport"      ,   &
835               "PW",1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
836            IF ( ln_ptrcomp ) THEN
837               CALL histdef( numptr, "sophtove", "Overturning Heat Transport"    ,   &
838                  "PW",1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
839            END IF
840            IF( ln_subbas ) THEN
841               CALL histdef( numptr, "sohtatl", "Heat Transport Atlantic"//TRIM(cl_comment),  &
842                  "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
843               CALL histdef( numptr, "sohtpac", "Heat Transport Pacific"//TRIM(cl_comment) ,  &
844                  "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
845               CALL histdef( numptr, "sohtind", "Heat Transport Indian"//TRIM(cl_comment)  ,  &
846                  "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
847               CALL histdef( numptr, "sohtipc", "Heat Transport Pacific+Indian"//TRIM(cl_comment), &
848                  "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
849            ENDIF
850
851
852            !  Salt transport
853
854            CALL histdef( numptr, "sopstadv", "Advective Salt Transport"      ,   &
855               "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
856            CALL histdef( numptr, "sopstldf", "Diffusive Salt Transport"      ,   &
857               "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
858            IF ( ln_ptrcomp ) THEN
859               CALL histdef( numptr, "sopstove", "Overturning Salt Transport"    ,   &
860                  "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
861            END IF
862#if defined key_diaeiv
863            ! Eddy induced velocity
864            CALL histdef( numptr, "zomsfeiv", "Bolus Meridional Stream-Function: global",   &
865               "Sv"      , 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout )
866            CALL histdef( numptr, "sophteiv", "Bolus Advective Heat Transport",   &
867               "PW"      , 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
868            CALL histdef( numptr, "sopsteiv", "Bolus Advective Salt Transport",   &
869               "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
870#endif
871            IF( ln_subbas ) THEN
872               CALL histdef( numptr, "sostatl", "Salt Transport Atlantic"//TRIM(cl_comment)      ,  &
873                  "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
874               CALL histdef( numptr, "sostpac", "Salt Transport Pacific"//TRIM(cl_comment)      ,   &
875                  "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
876               CALL histdef( numptr, "sostind", "Salt Transport Indian"//TRIM(cl_comment)      ,    &
877                  "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
878               CALL histdef( numptr, "sostipc", "Salt Transport Pacific+Indian"//TRIM(cl_comment),  &
879                  "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout )
880            ENDIF
881
882            CALL histend( numptr, snc4set )
883
884         END IF
885#if defined key_mpp_mpi
886      END IF
887#endif
888
889#if defined key_mpp_mpi
890      IF( MOD( itmod, nf_ptr ) == 0 .AND. l_znl_root ) THEN
891#else
892      IF( MOD( itmod, nf_ptr ) == 0  ) THEN
893#endif
894         niter = niter + 1
895
896!-$$         IF(lwp) THEN
897!-$$            WRITE(numout,*)
898!-$$            WRITE(numout,*) 'dia_ptr_wri : write Poleward Transports at time-step : kt = ', kt, &
899!-$$               & 'it = ', it, ' itmod = ', itmod, ' niter = ', niter
900!-$$            WRITE(numout,*) '~~~~~~~~~~'
901!-$$            WRITE(numout,*)
902!-$$         ENDIF
903
904         IF (ln_diaznl ) THEN
905            CALL histwrite( numptr, "zosrfglo", niter, surf_jk_glo , ndim, ndex )
906            CALL histwrite( numptr, "zotemglo", niter, tn_jk_glo  , ndim, ndex )
907            CALL histwrite( numptr, "zosalglo", niter, sn_jk_glo  , ndim, ndex )
908
909            IF (ln_subbas) THEN
910               CALL histwrite( numptr, "zosrfatl", niter, surf_jk_atl, ndim_atl, ndex_atl )
911               CALL histwrite( numptr, "zosrfpac", niter, surf_jk_pac, ndim_pac, ndex_pac )
912               CALL histwrite( numptr, "zosrfind", niter, surf_jk_ind, ndim_ind, ndex_ind )
913               CALL histwrite( numptr, "zosrfipc", niter, surf_jk_ipc, ndim_ipc, ndex_ipc )
914
915               CALL histwrite( numptr, "zotematl", niter, tn_jk_atl  , ndim_atl, ndex_atl )
916               CALL histwrite( numptr, "zosalatl", niter, sn_jk_atl  , ndim_atl, ndex_atl )
917               CALL histwrite( numptr, "zotempac", niter, tn_jk_pac  , ndim_pac, ndex_pac )
918               CALL histwrite( numptr, "zosalpac", niter, sn_jk_pac  , ndim_pac, ndex_pac )
919               CALL histwrite( numptr, "zotemind", niter, tn_jk_ind  , ndim_ind, ndex_ind )
920               CALL histwrite( numptr, "zosalind", niter, sn_jk_ind  , ndim_ind, ndex_ind )
921               CALL histwrite( numptr, "zotemipc", niter, tn_jk_ipc  , ndim_ipc, ndex_ipc )
922               CALL histwrite( numptr, "zosalipc", niter, sn_jk_ipc  , ndim_ipc, ndex_ipc )
923            END IF
924         ENDIF
925
926         ! overturning outputs:
927         CALL histwrite( numptr, "zomsfglo", niter, v_msf_glo, ndim, ndex )
928         IF( ln_subbas .AND. ln_diaznl ) THEN
929            CALL histwrite( numptr, "zomsfatl", niter, v_msf_atl , ndim_atl_30, ndex_atl_30 )
930            CALL histwrite( numptr, "zomsfpac", niter, v_msf_pac , ndim_pac_30, ndex_pac_30 )
931            CALL histwrite( numptr, "zomsfind", niter, v_msf_ind , ndim_ind_30, ndex_ind_30 )
932            CALL histwrite( numptr, "zomsfipc", niter, v_msf_ipc , ndim_ipc_30, ndex_ipc_30 )
933         ENDIF
934#if defined key_diaeiv
935         CALL histwrite( numptr, "zomsfeiv", niter, v_msf_eiv_glo, ndim  , ndex   )
936#endif
937
938
939         ! heat transport outputs:
940         IF( ln_subbas ) THEN
941            CALL histwrite( numptr, "sohtatl", niter, ht_atl  , ndim_h_atl_30, ndex_h_atl_30 )
942            CALL histwrite( numptr, "sohtpac", niter, ht_pac  , ndim_h_pac_30, ndex_h_pac_30 )
943            CALL histwrite( numptr, "sohtind", niter, ht_ind  , ndim_h_ind_30, ndex_h_ind_30 )
944            CALL histwrite( numptr, "sohtipc", niter, ht_ipc  , ndim_h_ipc_30, ndex_h_ipc_30 )
945            CALL histwrite( numptr, "sostatl", niter, st_atl  , ndim_h_atl_30, ndex_h_atl_30 )
946            CALL histwrite( numptr, "sostpac", niter, st_pac  , ndim_h_pac_30, ndex_h_pac_30 )
947            CALL histwrite( numptr, "sostind", niter, st_ind  , ndim_h_ind_30, ndex_h_ind_30 )
948            CALL histwrite( numptr, "sostipc", niter, st_ipc  , ndim_h_ipc_30, ndex_h_ipc_30 )
949         ENDIF
950
951         CALL histwrite( numptr, "sophtadv", niter, pht_adv     , ndim_h, ndex_h )
952         CALL histwrite( numptr, "sophtldf", niter, pht_ldf     , ndim_h, ndex_h )
953         CALL histwrite( numptr, "sopstadv", niter, pst_adv     , ndim_h, ndex_h )
954         CALL histwrite( numptr, "sopstldf", niter, pst_ldf     , ndim_h, ndex_h )
955         IF ( ln_ptrcomp ) THEN
956            CALL histwrite( numptr, "sopstove", niter, pst_ove_glo , ndim_h, ndex_h )
957            CALL histwrite( numptr, "sophtove", niter, pht_ove_glo , ndim_h, ndex_h )
958         ENDIF
959#if defined key_diaeiv
960         CALL histwrite( numptr, "sophteiv", niter, pht_eiv_glo  , ndim_h, ndex_h )
961         CALL histwrite( numptr, "sopsteiv", niter, pst_eiv_glo  , ndim_h, ndex_h )
962#endif
963         !
964      ENDIF
965      !
966   END SUBROUTINE dia_ptr_wri
967
968   !!======================================================================
969END MODULE diaptr
Note: See TracBrowser for help on using the repository browser.