MODULE diaptr !!====================================================================== !! *** MODULE diaptr *** !! Ocean physics: Computes meridonal transports and zonal means !!===================================================================== !! History : 9.0 ! 03-09 (C. Talandier, G. Madec) Original code !! 9.0 ! 06-01 (A. Biastoch) Allow sub-basins computation !! 9.0 ! 03-09 (O. Marti) Add fields !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! dia_ptr : Poleward Transport Diagnostics module !! dia_ptr_init : Initialization, namelist read !! dia_ptr_wri : Output of poleward fluxes !! ptr_vjk : "zonal" sum computation of a "meridional" flux array !! ptr_tjk : "zonal" mean computation of a tracer field !! ptr_vj : "zonal" and vertical sum computation of a "meridional" !! : flux array; Generic interface: ptr_vj_3d, ptr_vj_2d !!---------------------------------------------------------------------- USE oce ! ocean dynamics and active tracers USE dom_oce ! ocean space and time domain USE ldftra_oce ! ocean active tracers: lateral physics USE lib_mpp USE in_out_manager USE dianam USE phycst USE iom USE ioipsl USE daymod IMPLICIT NONE PRIVATE INTERFACE ptr_vj MODULE PROCEDURE ptr_vj_3d, ptr_vj_2d END INTERFACE PUBLIC dia_ptr_init ! call in opa module PUBLIC dia_ptr ! call in step module PUBLIC ptr_vj ! call by tra_ldf & tra_adv routines PUBLIC ptr_vjk ! call by tra_ldf & tra_adv routines !!! ** init namelist (namptr) LOGICAL , PUBLIC :: ln_diaptr = .FALSE. !: Poleward transport flag (T) or not (F) LOGICAL , PUBLIC :: ln_subbas = .FALSE. !: Atlantic/Pacific/Indian basins calculation LOGICAL , PUBLIC :: ln_diaznl = .FALSE. !: Add zonal means and meridional stream functions LOGICAL , PUBLIC :: ln_ptrcomp = .FALSE. !: Add decomposition : overturning (and gyre, soon ...) INTEGER , PUBLIC :: nf_ptr = 15 !: frequency of ptr computation INTEGER , PUBLIC :: nf_ptr_wri = 15 !: frequency of ptr outputs REAL(wp), DIMENSION(jpi,jpj), PUBLIC :: abasin, pbasin, ibasin, dbasin, sbasin !: Sub basin masks REAL(wp), PUBLIC, DIMENSION(jpj) :: pht_adv, pst_adv !: heat and salt poleward transport: advection REAL(wp), PUBLIC, DIMENSION(jpj) :: pht_ove_glo, pst_ove_glo, pht_ove_atl, pst_ove_atl, pht_ove_pac, pst_ove_pac, & & pht_ove_ind, pst_ove_ind, pht_ove_ipc, pst_ove_ipc !: heat and salt poleward transport: overturning REAL(wp), PUBLIC, DIMENSION(jpj) :: pht_ldf, pst_ldf !: heat and salt poleward transport: lateral diffusion #if defined key_diaeiv REAL(wp), PUBLIC, DIMENSION(jpj) :: pht_eiv_glo, pst_eiv_glo, pht_eiv_atl, pst_eiv_atl, pht_eiv_pac, pst_eiv_pac, & & pht_eiv_ind, pst_eiv_ind, pht_eiv_ipc, pst_eiv_ipc !: heat and salt poleward transport: bolus advection #endif REAL(wp), PUBLIC, DIMENSION(jpj) :: ht_glo,ht_atl,ht_ind,ht_pac,ht_ipc !: heat REAL(wp), PUBLIC, DIMENSION(jpj) :: st_glo,st_atl,st_ind,st_pac,st_ipc !: salt INTEGER :: niter INTEGER :: nidom_ptr REAL(wp), DIMENSION(jpj,jpk) :: tn_jk_glo, sn_jk_glo, & !: "zonal" mean temperature and salinity & tn_jk_atl, sn_jk_atl, & & tn_jk_pac, sn_jk_pac, & & tn_jk_ind, sn_jk_ind, & & tn_jk_ipc, sn_jk_ipc, & & v_msf_glo , & !: "meridional" Stream-Function & v_msf_atl , & & v_msf_pac , & & v_msf_ind , & & v_msf_ipc , & & surf_jk_glo , & !: Ocean "zonal" section surface & surf_jk_atl , & & surf_jk_pac , & & surf_jk_ind , & & surf_jk_ipc , & & surf_jk_r_glo , & !: inverse of the ocean "zonal" section surface & surf_jk_r_atl , & & surf_jk_r_pac , & & surf_jk_r_ind , & & surf_jk_r_ipc #if defined key_diaeiv REAL(wp), DIMENSION(jpj,jpk) :: v_msf_eiv_glo, v_msf_eiv_atl, v_msf_eiv_pac, & & v_msf_eiv_ind, v_msf_eiv_ipc !: bolus "meridional" Stream-Function #endif !! * Substitutions # include "domzgr_substitute.h90" # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! OPA 9.0 , LOCEAN-IPSL (2005) !! $Id$ !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS FUNCTION ptr_vj_3d( pva ) RESULT ( p_fval ) !!---------------------------------------------------------------------- !! *** ROUTINE ptr_vj_3d *** !! !! ** Purpose : "zonal" and vertical sum computation of a "meridional" !! flux array !! !! ** Method : - i-k sum of pva using the interior 2D vmask (vmask_i). !! pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v) !! !! ** Action : - p_fval: i-k-mean poleward flux of pva !!---------------------------------------------------------------------- REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk) :: pva ! mask flux array at V-point !! INTEGER :: ji, jj, jk ! dummy loop arguments INTEGER :: ijpj ! ??? REAL(wp), DIMENSION(jpj) :: p_fval ! function value !!-------------------------------------------------------------------- ! ijpj = jpj p_fval(:) = 0.e0 DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = fs_2, fs_jpim1 ! Vector opt. p_fval(jj) = p_fval(jj) + pva(ji,jj,jk) * tmask_i(ji,jj) END DO END DO END DO ! #if defined key_mpp_mpi CALL mpp_sum( p_fval, ijpj, ncomm_znl) !!bug I presume #endif ! END FUNCTION ptr_vj_3d FUNCTION ptr_vj_2d( pva ) RESULT ( p_fval ) !!---------------------------------------------------------------------- !! *** ROUTINE ptr_vj_2d *** !! !! ** Purpose : "zonal" and vertical sum computation of a "meridional" !! flux array !! !! ** Method : - i-k sum of pva using the interior 2D vmask (vmask_i). !! pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v) !! !! ** Action : - p_fval: i-k-mean poleward flux of pva !!---------------------------------------------------------------------- REAL(wp) , INTENT(in), DIMENSION(jpi,jpj) :: pva ! mask flux array at V-point !! INTEGER :: ji,jj ! dummy loop arguments INTEGER :: ijpj ! ??? REAL(wp), DIMENSION(jpj) :: p_fval ! function value !!-------------------------------------------------------------------- ! ijpj = jpj p_fval(:) = 0.e0 DO jj = 2, jpjm1 DO ji = nldi, nlei ! No vector optimisation here. Better use a mask ? p_fval(jj) = p_fval(jj) + pva(ji,jj) * tmask_i(ji,jj) END DO END DO ! #if defined key_mpp_mpi CALL mpp_sum( p_fval, ijpj, ncomm_znl ) !!bug I presume #endif ! END FUNCTION ptr_vj_2d FUNCTION ptr_vjk( pva, bmask ) RESULT ( p_fval ) !!---------------------------------------------------------------------- !! *** ROUTINE ptr_vjk *** !! !! ** Purpose : "zonal" sum computation of a "meridional" flux array !! !! ** Method : - i-sum of pva using the interior 2D vmask (vmask_i). !! pva is supposed to be a masked flux (i.e. * vmask*e1v*e3v) !! !! ** Action : - p_fval: i-k-mean poleward flux of pva !!---------------------------------------------------------------------- REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk) :: pva ! mask flux array at V-point REAL(wp) , INTENT(in), DIMENSION(jpi,jpj), OPTIONAL :: bmask ! Optional 2D basin mask !! INTEGER :: ji, jj, jk ! dummy loop arguments INTEGER , DIMENSION (1) :: ish INTEGER , DIMENSION (2) :: ish2 REAL(wp), DIMENSION(jpj*jpk) :: zwork ! temporary vector for mpp_sum REAL(wp), DIMENSION(jpj,jpk) :: p_fval ! return function value !!-------------------------------------------------------------------- ! p_fval(:,:) = 0.e0 ! IF (PRESENT (bmask)) THEN DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = nldi, nlei ! No vector optimisation here. Better use a mask ? p_fval(jj,jk) = p_fval(jj,jk) + pva(ji,jj,jk) * e1v(ji,jj) * fse3v(ji,jj,jk) & & * tmask_i(ji,jj) * bmask(ji,jj) END DO END DO END DO ELSE DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = nldi, nlei ! No vector optimisation here. Better use a mask ? p_fval(jj,jk) = p_fval(jj,jk) + pva(ji,jj,jk) * e1v(ji,jj) * fse3v(ji,jj,jk) & & * tmask_i(ji,jj) END DO END DO END DO END IF ! #if defined key_mpp_mpi ish(1) = jpj*jpk ; ish2(1) = jpj ; ish2(2) = jpk zwork(:)= RESHAPE( p_fval, ish ) CALL mpp_sum( zwork, jpj*jpk, ncomm_znl ) p_fval(:,:)= RESHAPE( zwork, ish2 ) #endif ! END FUNCTION ptr_vjk FUNCTION ptr_tjk( pta, bmask ) RESULT ( p_fval ) !!---------------------------------------------------------------------- !! *** ROUTINE ptr_tjk *** !! !! ** Purpose : "zonal" mean computation of a tracer field !! !! ** Method : - i-sum of mj(pta) using tmask !! multiplied by the inverse of the surface of the "zonal" ocean !! section !! !! ** Action : - p_fval: i-k-mean poleward flux of pta !!---------------------------------------------------------------------- REAL(wp) , INTENT(in), DIMENSION(jpi,jpj,jpk) :: pta ! tracer flux array at T-point REAL(wp) , INTENT(in), DIMENSION(jpi,jpj), OPTIONAL :: bmask ! Optional 2D basin mask !! INTEGER :: ji, jj, jk ! dummy loop arguments INTEGER, DIMENSION (1) :: ish INTEGER, DIMENSION (2) :: ish2 REAL(wp),DIMENSION(jpj*jpk) :: zwork ! temporary vector for mpp_sum REAL(wp),DIMENSION(jpj,jpk) :: p_fval ! return function value !!-------------------------------------------------------------------- ! p_fval(:,:) = 0.e0 IF (PRESENT (bmask)) THEN DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = nldi, nlei ! No vector optimisation here. Better use a mask ? p_fval(jj,jk) = p_fval(jj,jk) + pta(ji,jj,jk) & & * e1t(ji,jj) * fse3t(ji,jj,jk) & & * tmask_i(ji,jj) & & * bmask(ji,jj) END DO END DO END DO ELSE DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = nldi, nlei ! No vector optimisation here. Better use a mask ? p_fval(jj,jk) = p_fval(jj,jk) + pta(ji,jj,jk) & & * e1t(ji,jj) * fse3t(ji,jj,jk) & & * tmask_i(ji,jj) END DO END DO END DO END IF p_fval(:,:) = p_fval(:,:) * 0.5 #if defined key_mpp_mpi ish(1) = jpj*jpk ; ish2(1) = jpj ; ish2(2) = jpk zwork(:)= RESHAPE( p_fval, ish ) CALL mpp_sum( zwork, jpj*jpk, ncomm_znl ) p_fval(:,:)= RESHAPE(zwork,ish2) #endif ! END FUNCTION ptr_tjk SUBROUTINE dia_ptr( kt ) !!---------------------------------------------------------------------- !! *** ROUTINE dia_ptr *** !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time step index !! INTEGER :: jk, jj, ji ! dummy loop REAL(wp) :: zsverdrup, & ! conversion from m3/s to Sverdrup & zpwatt, & ! conversion from W to PW & zggram ! conversion from g to Pg REAL(wp), DIMENSION(jpi,jpj,jpk) :: vt, vs !!---------------------------------------------------------------------- IF( kt == nit000 .OR. MOD( kt, nf_ptr ) == 0 ) THEN IF ( MOD( kt, nf_ptr ) == 0 ) THEN zsverdrup = 1.e-6 zpwatt = 1.e-15 zggram = 1.e-6 IF ( ln_diaznl ) THEN ! "zonal" mean temperature and salinity at V-points tn_jk_glo(:,:) = ptr_tjk( tn(:,:,:) ) * surf_jk_r_glo(:,:) sn_jk_glo(:,:) = ptr_tjk( sn(:,:,:) ) * surf_jk_r_glo(:,:) IF (ln_subbas) THEN tn_jk_atl(:,:) = ptr_tjk( tn(:,:,:), abasin(:,:) ) * surf_jk_r_atl(:,:) sn_jk_atl(:,:) = ptr_tjk( sn(:,:,:), abasin(:,:) ) * surf_jk_r_atl(:,:) tn_jk_pac(:,:) = ptr_tjk( tn(:,:,:), pbasin(:,:) ) * surf_jk_r_pac(:,:) sn_jk_pac(:,:) = ptr_tjk( sn(:,:,:), pbasin(:,:) ) * surf_jk_r_pac(:,:) tn_jk_ind(:,:) = ptr_tjk( tn(:,:,:), ibasin(:,:) ) * surf_jk_r_ind(:,:) sn_jk_ind(:,:) = ptr_tjk( sn(:,:,:), ibasin(:,:) ) * surf_jk_r_ind(:,:) tn_jk_ipc(:,:) = ptr_tjk( tn(:,:,:), dbasin(:,:) ) * surf_jk_r_ipc(:,:) sn_jk_ipc(:,:) = ptr_tjk( sn(:,:,:), dbasin(:,:) ) * surf_jk_r_ipc(:,:) ENDIF ENDIF !-------------------------------------------------------- ! overturning calculation: ! horizontal integral and vertical dz #if defined key_diaeiv v_msf_glo(:,:) = ptr_vjk( vn(:,:,:)+v_eiv(:,:,:) ) IF( ln_subbas .AND. ln_diaznl ) THEN v_msf_atl(:,:) = ptr_vjk( vn (:,:,:)+v_eiv(:,:,:), abasin(:,:)*sbasin(:,:) ) v_msf_pac(:,:) = ptr_vjk( vn (:,:,:)+v_eiv(:,:,:), pbasin(:,:)*sbasin(:,:) ) v_msf_ind(:,:) = ptr_vjk( vn (:,:,:)+v_eiv(:,:,:), ibasin(:,:)*sbasin(:,:) ) v_msf_ipc(:,:) = ptr_vjk( vn (:,:,:)+v_eiv(:,:,:), dbasin(:,:)*sbasin(:,:) ) ENDIF #else v_msf_glo(:,:) = ptr_vjk( vn(:,:,:) ) IF( ln_subbas .AND. ln_diaznl ) THEN v_msf_atl(:,:) = ptr_vjk( vn (:,:,:), abasin(:,:)*sbasin(:,:) ) v_msf_pac(:,:) = ptr_vjk( vn (:,:,:), pbasin(:,:)*sbasin(:,:) ) v_msf_ind(:,:) = ptr_vjk( vn (:,:,:), ibasin(:,:)*sbasin(:,:) ) v_msf_ipc(:,:) = ptr_vjk( vn (:,:,:), dbasin(:,:)*sbasin(:,:) ) ENDIF #endif #if defined key_diaeiv v_msf_eiv_glo(:,:) = ptr_vjk( v_eiv(:,:,:) ) IF (ln_subbas ) THEN v_msf_eiv_atl(:,:) = ptr_vjk( v_eiv(:,:,:), abasin(:,:)*sbasin(:,:) ) v_msf_eiv_pac(:,:) = ptr_vjk( v_eiv(:,:,:), pbasin(:,:)*sbasin(:,:) ) v_msf_eiv_ind(:,:) = ptr_vjk( v_eiv(:,:,:), ibasin(:,:)*sbasin(:,:) ) v_msf_eiv_ipc(:,:) = ptr_vjk( v_eiv(:,:,:), dbasin(:,:)*sbasin(:,:) ) END IF #endif ! "Meridional" Stream-Function DO jk = 2,jpk v_msf_glo(:,jk) = v_msf_glo(:,jk-1) + v_msf_glo(:,jk) END DO v_msf_glo(:,:) = v_msf_glo(:,:) * zsverdrup #if defined key_diaeiv ! Bolus "Meridional" Stream-Function DO jk = 2,jpk v_msf_eiv_glo(:,jk) = v_msf_eiv_glo(:,jk-1) + v_msf_eiv_glo(:,jk) END DO v_msf_eiv_glo(:,:) = v_msf_eiv_glo(:,:) * zsverdrup IF ( ln_subbas ) THEN DO jk = 2,jpk v_msf_eiv_atl(:,jk) = v_msf_eiv_atl(:,jk-1) + v_msf_eiv_atl(:,jk) v_msf_eiv_pac(:,jk) = v_msf_eiv_pac(:,jk-1) + v_msf_eiv_pac(:,jk) v_msf_eiv_ind(:,jk) = v_msf_eiv_ind(:,jk-1) + v_msf_eiv_ind(:,jk) v_msf_eiv_ipc(:,jk) = v_msf_eiv_ipc(:,jk-1) + v_msf_eiv_ipc(:,jk) END DO ENDIF #endif ! IF( ln_subbas .AND. ln_diaznl ) THEN DO jk = 2,jpk v_msf_atl(:,jk) = v_msf_atl(:,jk-1) + v_msf_atl(:,jk) v_msf_pac(:,jk) = v_msf_pac(:,jk-1) + v_msf_pac(:,jk) v_msf_ind(:,jk) = v_msf_ind(:,jk-1) + v_msf_ind(:,jk) v_msf_ipc(:,jk) = v_msf_ipc(:,jk-1) + v_msf_ipc(:,jk) END DO v_msf_atl(:,:) = v_msf_atl(:,:) * zsverdrup v_msf_pac(:,:) = v_msf_pac(:,:) * zsverdrup v_msf_ind(:,:) = v_msf_ind(:,:) * zsverdrup v_msf_ipc(:,:) = v_msf_ipc(:,:) * zsverdrup ENDIF ! Transports ! T times V on T points (include bolus velocities) #if defined key_diaeiv DO jj = 1, jpj DO ji = 1, jpi vt(ji,jj,:) = tn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) + u_eiv(ji,jj,:) + u_eiv(ji,jj-1,:) )*0.5 vs(ji,jj,:) = sn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) + v_eiv(ji,jj,:) + v_eiv(ji,jj-1,:) )*0.5 END DO END DO #else DO jj = 1, jpj DO ji = 1, jpi vt(ji,jj,:) = tn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) )*0.5 vs(ji,jj,:) = sn(ji,jj,:) * ( vn(ji,jj,:) + vn(ji,jj-1,:) )*0.5 END DO END DO #endif ht_glo(:) = SUM( ptr_vjk( vt(:,:,:)), 2 ) st_glo(:) = SUM( ptr_vjk( vs(:,:,:)), 2 ) IF ( ln_subbas ) THEN ht_atl(:) = SUM( ptr_vjk( vt (:,:,:), abasin(:,:)*sbasin(:,:)), 2 ) ht_pac(:) = SUM( ptr_vjk( vt (:,:,:), pbasin(:,:)*sbasin(:,:)), 2 ) ht_ind(:) = SUM( ptr_vjk( vt (:,:,:), ibasin(:,:)*sbasin(:,:)), 2 ) ht_ipc(:) = SUM( ptr_vjk( vt (:,:,:), dbasin(:,:)*sbasin(:,:)), 2 ) st_atl(:) = SUM( ptr_vjk( vs (:,:,:), abasin(:,:)*sbasin(:,:)), 2 ) st_pac(:) = SUM( ptr_vjk( vs (:,:,:), pbasin(:,:)*sbasin(:,:)), 2 ) st_ind(:) = SUM( ptr_vjk( vs (:,:,:), ibasin(:,:)*sbasin(:,:)), 2 ) st_ipc(:) = SUM( ptr_vjk( vs (:,:,:), dbasin(:,:)*sbasin(:,:)), 2 ) ENDIF ! poleward tracer transports: ! overturning components: IF ( ln_ptrcomp ) THEN pht_ove_glo(:) = SUM( v_msf_glo(:,:) * tn_jk_glo(:,:), 2 ) ! SUM over jk pst_ove_glo(:) = SUM( v_msf_glo(:,:) * sn_jk_glo(:,:), 2 ) IF ( ln_subbas ) THEN pht_ove_atl(:) = SUM( v_msf_atl(:,:) * tn_jk_atl(:,:), 2 ) ! SUM over jk pst_ove_atl(:) = SUM( v_msf_atl(:,:) * sn_jk_atl(:,:), 2 ) pht_ove_pac(:) = SUM( v_msf_pac(:,:) * tn_jk_pac(:,:), 2 ) ! SUM over jk pst_ove_pac(:) = SUM( v_msf_pac(:,:) * sn_jk_pac(:,:), 2 ) pht_ove_ind(:) = SUM( v_msf_ind(:,:) * tn_jk_ind(:,:), 2 ) ! SUM over jk pst_ove_ind(:) = SUM( v_msf_ind(:,:) * sn_jk_ind(:,:), 2 ) pht_ove_ipc(:) = SUM( v_msf_ipc(:,:) * tn_jk_ipc(:,:), 2 ) ! SUM over jk pst_ove_ipc(:) = SUM( v_msf_ipc(:,:) * sn_jk_ipc(:,:), 2 ) END IF END IF ! Bolus component #if defined key_diaeiv pht_eiv_glo(:) = SUM( v_msf_eiv_glo(:,:) * tn_jk_glo(:,:), 2 ) ! SUM over jk pst_eiv_glo(:) = SUM( v_msf_eiv_glo(:,:) * sn_jk_glo(:,:), 2 ) ! SUM over jk IF ( ln_subbas ) THEN pht_eiv_atl(:) = SUM( v_msf_eiv_glo(:,:) * tn_jk_atl(:,:), 2 ) ! SUM over jk pst_eiv_atl(:) = SUM( v_msf_eiv_glo(:,:) * sn_jk_atl(:,:), 2 ) ! SUM over jk pht_eiv_pac(:) = SUM( v_msf_eiv_pac(:,:) * tn_jk_pac(:,:), 2 ) ! SUM over jk pst_eiv_pac(:) = SUM( v_msf_eiv_pac(:,:) * sn_jk_pac(:,:), 2 ) ! SUM over jk pht_eiv_ind(:) = SUM( v_msf_eiv_ind(:,:) * tn_jk_ind(:,:), 2 ) ! SUM over jk pst_eiv_ind(:) = SUM( v_msf_eiv_ind(:,:) * sn_jk_ind(:,:), 2 ) ! SUM over jk pht_eiv_ipc(:) = SUM( v_msf_eiv_ipc(:,:) * tn_jk_ipc(:,:), 2 ) ! SUM over jk pst_eiv_ipc(:) = SUM( v_msf_eiv_ipc(:,:) * sn_jk_ipc(:,:), 2 ) ! SUM over jk ENDIF #endif ! conversion in PW and G g zpwatt = zpwatt * rau0 * rcp pht_adv(:) = pht_adv(:) * zpwatt pht_ldf(:) = pht_ldf(:) * zpwatt pst_adv(:) = pst_adv(:) * zggram pst_ldf(:) = pst_ldf(:) * zggram IF ( ln_ptrcomp ) THEN pht_ove_glo(:) = pht_ove_glo(:) * zpwatt pst_ove_glo(:) = pst_ove_glo(:) * zggram END IF #if defined key_diaeiv pht_eiv_glo(:) = pht_eiv_glo(:) * zpwatt pst_eiv_glo(:) = pst_eiv_glo(:) * zggram #endif IF( ln_subbas ) THEN ht_atl(:) = ht_atl(:) * zpwatt ht_pac(:) = ht_pac(:) * zpwatt ht_ind(:) = ht_ind(:) * zpwatt ht_ipc(:) = ht_ipc(:) * zpwatt st_atl(:) = st_atl(:) * zggram st_pac(:) = st_pac(:) * zggram st_ind(:) = st_ind(:) * zggram st_ipc(:) = st_ipc(:) * zggram ENDIF ENDIF ! outputs CALL dia_ptr_wri( kt ) ENDIF ! Close the file IF( kt == nitend ) CALL histclo( numptr ) ! END SUBROUTINE dia_ptr SUBROUTINE dia_ptr_init !!---------------------------------------------------------------------- !! *** ROUTINE dia_ptr_init *** !! !! ** Purpose : Initialization, namelist read !!---------------------------------------------------------------------- NAMELIST/namptr/ ln_diaptr, ln_diaznl, ln_subbas, ln_ptrcomp, nf_ptr, nf_ptr_wri INTEGER :: inum ! temporary logical unit #if defined key_mpp_mpi INTEGER, DIMENSION (1) :: iglo, iloc, iabsf, iabsl, ihals, ihale, idid #endif !!---------------------------------------------------------------------- ! Read Namelist namptr : poleward transport parameters REWIND ( numnam ) READ ( numnam, namptr ) ! Control print IF(lwp) THEN WRITE(numout,*) WRITE(numout,*) 'dia_ptr_init : poleward transport and msf initialization' WRITE(numout,*) '~~~~~~~~~~~~' WRITE(numout,*) ' Namelist namptr : set ptr parameters' WRITE(numout,*) ' Switch for ptr diagnostic (T) or not (F) ln_diaptr = ', ln_diaptr WRITE(numout,*) ' Atla/Paci/Ind basins computation ln_subbas = ', ln_subbas WRITE(numout,*) ' Frequency of computation nf_ptr = ', nf_ptr WRITE(numout,*) ' Frequency of outputs nf_ptr_wri = ', nf_ptr_wri ENDIF ! ! Define MPI communicator for zonal sum ! IF( lk_mpp ) THEN CALL mpp_ini_znl ENDIF IF( ln_subbas ) THEN ! load sub-basin mask CALL iom_open( 'subbasins', inum ) CALL iom_get( inum, jpdom_data, 'atlmsk', abasin ) ! Atlantic basin CALL iom_get( inum, jpdom_data, 'pacmsk', pbasin ) ! Pacific basin CALL iom_get( inum, jpdom_data, 'indmsk', ibasin ) ! Indian basin CALL iom_close( inum ) dbasin(:,:) = MAX ( pbasin(:,:), ibasin(:,:) ) sbasin(:,:) = tmask (:,:,1) WHERE ( gphit (:,:) < -30.e0) sbasin(:,:) = 0.e0 ENDIF ! inverse of the ocean "zonal" v-point section surf_jk_glo(:,:) = ptr_tjk( tmask(:,:,:) ) surf_jk_r_glo(:,:) = 0.e0 WHERE( surf_jk_glo(:,:) /= 0.e0 ) surf_jk_r_glo(:,:) = 1.e0 / surf_jk_glo(:,:) IF (ln_subbas) THEN surf_jk_atl(:,:) = ptr_tjk( tmask (:,:,:), abasin(:,:) ) surf_jk_r_atl(:,:) = 0.e0 WHERE( surf_jk_atl(:,:) /= 0.e0 ) surf_jk_r_atl(:,:) = 1.e0 / surf_jk_atl(:,:) ! surf_jk_pac(:,:) = ptr_tjk( tmask (:,:,:), pbasin(:,:) ) surf_jk_r_pac(:,:) = 0.e0 WHERE( surf_jk_pac(:,:) /= 0.e0 ) surf_jk_r_pac(:,:) = 1.e0 / surf_jk_pac(:,:) ! surf_jk_ind(:,:) = ptr_tjk( tmask (:,:,:), ibasin(:,:) ) surf_jk_r_ind(:,:) = 0.e0 WHERE( surf_jk_ind(:,:) /= 0.e0 ) surf_jk_r_ind(:,:) = 1.e0 / surf_jk_ind(:,:) ! surf_jk_ipc(:,:) = ptr_tjk( tmask (:,:,:), dbasin(:,:) ) surf_jk_r_ipc(:,:) = 0.e0 WHERE( surf_jk_ipc(:,:) /= 0.e0 ) surf_jk_r_ipc(:,:) = 1.e0 / surf_jk_ipc(:,:) END IF !!---------------------------------------------------------------------- #if defined key_mpp_mpi iglo (1) = jpjglo iloc (1) = nlcj iabsf(1) = njmppt(narea) iabsl(:) = iabsf(:) + iloc(:) - 1 ihals(1) = nldj - 1 ihale(1) = nlcj - nlej idid (1) = 2 !-$$ IF(lwp) THEN !-$$ WRITE(numout,*) !-$$ WRITE(numout,*) 'dia_ptr_init : iloc = ', iloc !-$$ WRITE(numout,*) '~~~~~~~~~~~~ iabsf = ', iabsf !-$$ WRITE(numout,*) ' ihals = ', ihals !-$$ WRITE(numout,*) ' ihale = ', ihale !-$$ ENDIF CALL flio_dom_set ( jpnj, nproc/jpni, idid, iglo, iloc, iabsf, iabsl, ihals, ihale, 'BOX', nidom_ptr) #else nidom_ptr = FLIO_DOM_NONE #endif END SUBROUTINE dia_ptr_init SUBROUTINE dia_ptr_wri( kt ) !!--------------------------------------------------------------------- !! *** ROUTINE dia_ptr_wri *** !! !! ** Purpose : output of poleward fluxes !! !! ** Method : NetCDF file !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time-step index !! INTEGER, SAVE :: nhoridz, ndepidzt, ndepidzw INTEGER, SAVE :: ndim , ndim_atl , ndim_pac , ndim_ind , ndim_ipc INTEGER, SAVE :: ndim_atl_30 , ndim_pac_30 , ndim_ind_30 , ndim_ipc_30 INTEGER, SAVE :: ndim_h, ndim_h_atl_30, ndim_h_pac_30, ndim_h_ind_30, ndim_h_ipc_30 INTEGER, SAVE, DIMENSION (jpj*jpk) :: ndex , ndex_atl , ndex_pac , ndex_ind , ndex_ipc INTEGER, SAVE, DIMENSION (jpj*jpk) :: ndex_atl_30 , ndex_pac_30 , ndex_ind_30 , ndex_ipc_30 INTEGER, SAVE, DIMENSION (jpj) :: ndex_h, ndex_h_atl_30, ndex_h_pac_30, ndex_h_ind_30, ndex_h_ipc_30 CHARACTER (len=40) :: clhstnam, clop, clop_once, cl_comment ! temporary names INTEGER :: iline, it, itmod, ji, jj, jk ! REAL(wp) :: zsto, zout, zdt, zjulian ! temporary scalars REAL(wp), DIMENSION(jpj) :: zphi, zfoo REAL(wp), DIMENSION(jpj,jpk) :: z_1 !!---------------------------------------------------------------------- ! define time axis it = kt / nf_ptr itmod = kt - nit000 + 1 !-$$ IF(lwp) THEN !-$$ WRITE(numout,*) !-$$ WRITE(numout,*) 'dia_ptr_wri : kt = ', kt, 'it = ', it, ' itmod = ', itmod, ' niter = ', niter !-$$ WRITE(numout,*) '~~~~~~~~~~~~' !-$$ ENDIF ! Initialization ! -------------- IF( kt == nit000 ) THEN niter = (nit000 - 1) / nf_ptr !-$$ IF(lwp) THEN !-$$ WRITE(numout,*) !-$$ WRITE(numout,*) 'dia_ptr_wri : poleward transport and msf writing: initialization , niter = ', niter !-$$ WRITE(numout,*) '~~~~~~~~~~~~' !-$$ ENDIF zdt = rdt IF( nacc == 1 ) zdt = rdtmin ! Reference latitude ! ------------------ ! ! ======================= IF( cp_cfg == "orca" ) THEN ! ORCA configurations ! ! ======================= IF( jp_cfg == 05 ) iline = 192 ! i-line that passes near the North Pole IF( jp_cfg == 025 ) iline = 384 ! i-line that passes near the North Pole IF( jp_cfg == 1 ) iline = 96 ! i-line that passes near the North Pole IF( jp_cfg == 2 ) iline = 48 ! i-line that passes near the North Pole IF( jp_cfg == 4 ) iline = 24 ! i-line that passes near the North Pole zphi(:) = 0.e0 DO ji = mi0(iline), mi1(iline) zphi(:) = gphiv(ji,:) ! if iline is in the local domain ! Correct highest latitude for some configurations - will work if domain is parallelized in J ? IF( jp_cfg == 05 ) THEN DO jj = mj0(jpjdta), mj1(jpjdta) zphi( jj ) = zphi(mj0(jpjdta-1)) + (zphi(mj0(jpjdta-1))-zphi(mj0(jpjdta-2)))/2. zphi( jj ) = MIN( zphi(jj), 90.) END DO END IF IF( jp_cfg == 1 .OR. jp_cfg == 2 .OR. jp_cfg == 4 ) THEN DO jj = mj0(jpjdta-1), mj1(jpjdta-1) zphi( jj ) = 88.5e0 END DO DO jj = mj0(jpjdta ), mj1(jpjdta ) zphi( jj ) = 89.5e0 END DO END IF END DO ! provide the correct zphi to all local domains #if defined key_mpp_mpi CALL mpp_sum( zphi, jpj, ncomm_znl ) #endif ! ! ======================= ELSE ! OTHER configurations zjulian = zjulian - adatrj ! set calendar origin to the beginning of the experiment ! ! ======================= zphi(:) = gphiv(1,:) ! assume lat/lon coordinate, select the first i-line ! ENDIF ! ! Work only on westmost processor (will not work if mppini2 is used) #if defined key_mpp_mpi IF ( l_znl_root ) THEN #endif ! ! OPEN netcdf file ! ---------------- ! Define frequency of output and means zsto = nf_ptr * zdt IF( ln_mskland ) THEN ! put 1.e+20 on land (very expensive!!) clop = "ave(only(x))" clop_once = "once(only(x))" ELSE ! no use of the mask value (require less cpu time) clop = "ave(x)" clop_once = "once" ENDIF zout = nf_ptr_wri * zdt zfoo(:) = 0.e0 ! Compute julian date from starting date of the run CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian ) zjulian = zjulian - adatrj ! set calendar origin to the beginning of the experiment CALL dia_nam( clhstnam, nf_ptr_wri, 'diaptr' ) IF(lwp)WRITE( numout,*)" Name of diaptr NETCDF file : ", clhstnam ! Horizontal grid : zphi() CALL histbeg(clhstnam, 1, zfoo, jpj, zphi, & 1, 1, 1, jpj, niter, zjulian, zdt*nf_ptr, nhoridz, numptr, domain_id=nidom_ptr) ! Vertical grids : gdept_0, gdepw_0 CALL histvert( numptr, "deptht", "Vertical T levels", & "m", jpk, gdept_0, ndepidzt, "down" ) CALL histvert( numptr, "depthw", "Vertical W levels", & "m", jpk, gdepw_0, ndepidzw, "down" ) ! CALL wheneq ( jpj*jpk, MIN(surf_jk_glo(:,:), 1.e0), 1, 1., ndex , ndim ) ! Lat-Depth CALL wheneq ( jpj , MIN(surf_jk_glo(:,1), 1.e0), 1, 1., ndex_h, ndim_h ) ! Lat IF (ln_subbas) THEN z_1 (:,1) = 1.0e0 WHERE ( gphit (jpi/2,:) .LT. -30 ) z_1 (:,1) = 0.e0 DO jk = 2, jpk z_1 (:,jk) = z_1 (:,1) END DO CALL wheneq ( jpj*jpk, MIN(surf_jk_atl(:,:) , 1.e0), 1, 1., ndex_atl , ndim_atl ) ! Lat-Depth CALL wheneq ( jpj*jpk, MIN(surf_jk_atl(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_atl_30 , ndim_atl_30 ) ! Lat-Depth CALL wheneq ( jpj , MIN(surf_jk_atl(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_atl_30, ndim_h_atl_30 ) ! Lat CALL wheneq ( jpj*jpk, MIN(surf_jk_pac(:,:) , 1.e0), 1, 1., ndex_pac , ndim_pac ) ! Lat-Depth CALL wheneq ( jpj*jpk, MIN(surf_jk_pac(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_pac_30 , ndim_pac_30 ) ! Lat-Depth CALL wheneq ( jpj , MIN(surf_jk_pac(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_pac_30, ndim_h_pac_30 ) ! Lat CALL wheneq ( jpj*jpk, MIN(surf_jk_ind(:,:) , 1.e0), 1, 1., ndex_ind , ndim_ind ) ! Lat-Depth CALL wheneq ( jpj*jpk, MIN(surf_jk_ind(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_ind_30 , ndim_ind_30 ) ! Lat-Depth CALL wheneq ( jpj , MIN(surf_jk_ind(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_ind_30, ndim_h_ind_30 ) ! Lat CALL wheneq ( jpj*jpk, MIN(surf_jk_ipc(:,:) , 1.e0), 1, 1., ndex_ipc , ndim_ipc ) ! Lat-Depth CALL wheneq ( jpj*jpk, MIN(surf_jk_ipc(:,:)*z_1(:,:), 1.e0), 1, 1., ndex_ipc_30 , ndim_ipc_30 ) ! Lat-Depth CALL wheneq ( jpj , MIN(surf_jk_ipc(:,1)*z_1(:,1), 1.e0), 1, 1., ndex_h_ipc_30, ndim_h_ipc_30 ) ! Lat ENDIF ! #if defined key_diaeiv cl_comment = ' (Bolus part included)' #else cl_comment = ' ' #endif ! Zonal mean T and S IF ( ln_diaznl ) THEN CALL histdef( numptr, "zotemglo", "Zonal Mean Temperature","C" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosalglo", "Zonal Mean Salinity","PSU" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosrfglo", "Zonal Mean Surface","m^2" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout ) IF (ln_subbas) THEN CALL histdef( numptr, "zotematl", "Zonal Mean Temperature: Atlantic","C" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosalatl", "Zonal Mean Salinity: Atlantic","PSU" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosrfatl", "Zonal Mean Surface: Atlantic","m^2" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout ) CALL histdef( numptr, "zotempac", "Zonal Mean Temperature: Pacific","C" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosalpac", "Zonal Mean Salinity: Pacific","PSU" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosrfpac", "Zonal Mean Surface: Pacific","m^2" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout ) CALL histdef( numptr, "zotemind", "Zonal Mean Temperature: Indian","C" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosalind", "Zonal Mean Salinity: Indian","PSU" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosrfind", "Zonal Mean Surface: Indian","m^2" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout ) CALL histdef( numptr, "zotemipc", "Zonal Mean Temperature: Pacific+Indian","C" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosalipc", "Zonal Mean Salinity: Pacific+Indian","PSU" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop, zsto, zout ) CALL histdef( numptr, "zosrfipc", "Zonal Mean Surface: Pacific+Indian","m^2" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzt, 32, clop_once, zsto, zout ) ENDIF ENDIF ! Meridional Stream-Function (Eulerian and Bolus) CALL histdef( numptr, "zomsfglo", "Meridional Stream-Function: Global"//TRIM(cl_comment),"Sv" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) IF( ln_subbas .AND. ln_diaznl ) THEN CALL histdef( numptr, "zomsfatl", "Meridional Stream-Function: Atlantic"//TRIM(cl_comment),"Sv" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) CALL histdef( numptr, "zomsfpac", "Meridional Stream-Function: Pacific"//TRIM(cl_comment),"Sv" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) CALL histdef( numptr, "zomsfind", "Meridional Stream-Function: Indian"//TRIM(cl_comment),"Sv" , & 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) CALL histdef( numptr, "zomsfipc", "Meridional Stream-Function: Indo-Pacific"//TRIM(cl_comment),"Sv" ,& 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) ENDIF ! Heat transport CALL histdef( numptr, "sophtadv", "Advective Heat Transport" , & "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sophtldf", "Diffusive Heat Transport" , & "PW",1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) IF ( ln_ptrcomp ) THEN CALL histdef( numptr, "sophtove", "Overturning Heat Transport" , & "PW",1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) END IF IF( ln_subbas ) THEN CALL histdef( numptr, "sohtatl", "Heat Transport Atlantic"//TRIM(cl_comment), & "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sohtpac", "Heat Transport Pacific"//TRIM(cl_comment) , & "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sohtind", "Heat Transport Indian"//TRIM(cl_comment) , & "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sohtipc", "Heat Transport Pacific+Indian"//TRIM(cl_comment), & "PW", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) ENDIF ! Salt transport CALL histdef( numptr, "sopstadv", "Advective Salt Transport" , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sopstldf", "Diffusive Salt Transport" , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) IF ( ln_ptrcomp ) THEN CALL histdef( numptr, "sopstove", "Overturning Salt Transport" , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) END IF #if defined key_diaeiv ! Eddy induced velocity CALL histdef( numptr, "zomsfeiv", "Bolus Meridional Stream-Function: global", & "Sv" , 1, jpj, nhoridz, jpk, 1, jpk, ndepidzw, 32, clop, zsto, zout ) CALL histdef( numptr, "sophteiv", "Bolus Advective Heat Transport", & "PW" , 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sopsteiv", "Bolus Advective Salt Transport", & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) #endif IF( ln_subbas ) THEN CALL histdef( numptr, "sostatl", "Salt Transport Atlantic"//TRIM(cl_comment) , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sostpac", "Salt Transport Pacific"//TRIM(cl_comment) , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sostind", "Salt Transport Indian"//TRIM(cl_comment) , & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) CALL histdef( numptr, "sostipc", "Salt Transport Pacific+Indian"//TRIM(cl_comment), & "Giga g/s", 1, jpj, nhoridz, 1, 1, 1, -99, 32, clop, zsto, zout ) ENDIF CALL histend( numptr ) END IF #if defined key_mpp_mpi END IF #endif #if defined key_mpp_mpi IF( MOD( itmod, nf_ptr ) == 0 .AND. l_znl_root ) THEN #else IF( MOD( itmod, nf_ptr ) == 0 ) THEN #endif niter = niter + 1 !-$$ IF(lwp) THEN !-$$ WRITE(numout,*) !-$$ WRITE(numout,*) 'dia_ptr_wri : write Poleward Transports at time-step : kt = ', kt, & !-$$ & 'it = ', it, ' itmod = ', itmod, ' niter = ', niter !-$$ WRITE(numout,*) '~~~~~~~~~~' !-$$ WRITE(numout,*) !-$$ ENDIF IF (ln_diaznl ) THEN CALL histwrite( numptr, "zosrfglo", niter, surf_jk_glo , ndim, ndex ) CALL histwrite( numptr, "zotemglo", niter, tn_jk_glo , ndim, ndex ) CALL histwrite( numptr, "zosalglo", niter, sn_jk_glo , ndim, ndex ) IF (ln_subbas) THEN CALL histwrite( numptr, "zosrfatl", niter, surf_jk_atl, ndim_atl, ndex_atl ) CALL histwrite( numptr, "zosrfpac", niter, surf_jk_pac, ndim_pac, ndex_pac ) CALL histwrite( numptr, "zosrfind", niter, surf_jk_ind, ndim_ind, ndex_ind ) CALL histwrite( numptr, "zosrfipc", niter, surf_jk_ipc, ndim_ipc, ndex_ipc ) CALL histwrite( numptr, "zotematl", niter, tn_jk_atl , ndim_atl, ndex_atl ) CALL histwrite( numptr, "zosalatl", niter, sn_jk_atl , ndim_atl, ndex_atl ) CALL histwrite( numptr, "zotempac", niter, tn_jk_pac , ndim_pac, ndex_pac ) CALL histwrite( numptr, "zosalpac", niter, sn_jk_pac , ndim_pac, ndex_pac ) CALL histwrite( numptr, "zotemind", niter, tn_jk_ind , ndim_ind, ndex_ind ) CALL histwrite( numptr, "zosalind", niter, sn_jk_ind , ndim_ind, ndex_ind ) CALL histwrite( numptr, "zotemipc", niter, tn_jk_ipc , ndim_ipc, ndex_ipc ) CALL histwrite( numptr, "zosalipc", niter, sn_jk_ipc , ndim_ipc, ndex_ipc ) END IF ENDIF ! overturning outputs: CALL histwrite( numptr, "zomsfglo", niter, v_msf_glo, ndim, ndex ) IF( ln_subbas .AND. ln_diaznl ) THEN CALL histwrite( numptr, "zomsfatl", niter, v_msf_atl , ndim_atl_30, ndex_atl_30 ) CALL histwrite( numptr, "zomsfpac", niter, v_msf_pac , ndim_pac_30, ndex_pac_30 ) CALL histwrite( numptr, "zomsfind", niter, v_msf_ind , ndim_ind_30, ndex_ind_30 ) CALL histwrite( numptr, "zomsfipc", niter, v_msf_ipc , ndim_ipc_30, ndex_ipc_30 ) ENDIF #if defined key_diaeiv CALL histwrite( numptr, "zomsfeiv", niter, v_msf_eiv_glo, ndim , ndex ) #endif ! heat transport outputs: IF( ln_subbas ) THEN CALL histwrite( numptr, "sohtatl", niter, ht_atl , ndim_h_atl_30, ndex_h_atl_30 ) CALL histwrite( numptr, "sohtpac", niter, ht_pac , ndim_h_pac_30, ndex_h_pac_30 ) CALL histwrite( numptr, "sohtind", niter, ht_ind , ndim_h_ind_30, ndex_h_ind_30 ) CALL histwrite( numptr, "sohtipc", niter, ht_ipc , ndim_h_ipc_30, ndex_h_ipc_30 ) CALL histwrite( numptr, "sostatl", niter, st_atl , ndim_h_atl_30, ndex_h_atl_30 ) CALL histwrite( numptr, "sostpac", niter, st_pac , ndim_h_pac_30, ndex_h_pac_30 ) CALL histwrite( numptr, "sostind", niter, st_ind , ndim_h_ind_30, ndex_h_ind_30 ) CALL histwrite( numptr, "sostipc", niter, st_ipc , ndim_h_ipc_30, ndex_h_ipc_30 ) ENDIF CALL histwrite( numptr, "sophtadv", niter, pht_adv , ndim_h, ndex_h ) CALL histwrite( numptr, "sophtldf", niter, pht_ldf , ndim_h, ndex_h ) CALL histwrite( numptr, "sopstadv", niter, pst_adv , ndim_h, ndex_h ) CALL histwrite( numptr, "sopstldf", niter, pst_ldf , ndim_h, ndex_h ) IF ( ln_ptrcomp ) THEN CALL histwrite( numptr, "sopstove", niter, pst_ove_glo , ndim_h, ndex_h ) CALL histwrite( numptr, "sophtove", niter, pht_ove_glo , ndim_h, ndex_h ) ENDIF #if defined key_diaeiv CALL histwrite( numptr, "sophteiv", niter, pht_eiv_glo , ndim_h, ndex_h ) CALL histwrite( numptr, "sopsteiv", niter, pst_eiv_glo , ndim_h, ndex_h ) #endif ! ENDIF ! END SUBROUTINE dia_ptr_wri !!====================================================================== END MODULE diaptr