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

source: trunk/NEMO/OPA_SRC/DIA/diaptr.F90 @ 1413

Last change on this file since 1413 was 1413, checked in by flavoni, 15 years ago

correct bug for Meridionnal Stream function in diaptr.F90 when using mpi, see ticket #421 and #361

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