Changeset 4990 for trunk/NEMOGCM/NEMO/OPA_SRC/TRD
- Timestamp:
- 2014-12-15T17:42:49+01:00 (9 years ago)
- Location:
- trunk/NEMOGCM/NEMO/OPA_SRC/TRD
- Files:
-
- 8 deleted
- 3 edited
- 10 copied
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMOGCM/NEMO/OPA_SRC/TRD/trdtra.F90
r3632 r4990 2 2 !!====================================================================== 3 3 !! *** MODULE trdtra *** 4 !! Ocean diagnostics: ocean tracers trends 4 !! Ocean diagnostics: ocean tracers trends pre-processing 5 5 !!===================================================================== 6 !! History : 1.0 ! 2004-08 (C. Talandier) Original code 7 !! 2.0 ! 2005-04 (C. Deltel) Add Asselin trend in the ML budget 8 !! 3.3 ! 2010-06 (C. Ethe) merge TRA-TRC 9 !!---------------------------------------------------------------------- 10 #if defined key_trdtra || defined key_trdtrc || defined key_trdmld || defined key_trdmld_trc 11 !!---------------------------------------------------------------------- 12 !! trd_tra : Call the trend to be computed 13 !!---------------------------------------------------------------------- 14 USE dom_oce ! ocean domain 15 USE trdmod_oce ! ocean active mixed layer tracers trends 16 USE trdmod ! ocean active mixed layer tracers trends 17 USE trdmod_trc ! ocean passive mixed layer tracers trends 18 USE in_out_manager ! I/O manager 19 USE lib_mpp ! MPP library 20 USE wrk_nemo ! Memory allocation 21 6 !! History : 3.3 ! 2010-06 (C. Ethe) creation for the TRA/TRC merge 7 !! 3.5 ! 2012-02 (G. Madec) update the comments 8 !!---------------------------------------------------------------------- 9 10 !!---------------------------------------------------------------------- 11 !! trd_tra : pre-process the tracer trends 12 !! trd_tra_adv : transform a div(U.T) trend into a U.grad(T) trend 13 !! trd_tra_mng : tracer trend manager: dispatch to the diagnostic modules 14 !! trd_tra_iom : output 3D tracer trends using IOM 15 !!---------------------------------------------------------------------- 16 USE oce ! ocean dynamics and tracers variables 17 USE dom_oce ! ocean domain 18 USE sbc_oce ! surface boundary condition: ocean 19 USE zdf_oce ! ocean vertical physics 20 USE trd_oce ! trends: ocean variables 21 USE trdtrc ! ocean passive mixed layer tracers trends 22 USE trdglo ! trends: global domain averaged 23 USE trdpen ! trends: Potential ENergy 24 USE trdmxl ! ocean active mixed layer tracers trends 25 USE ldftra_oce ! ocean active tracers lateral physics 26 USE zdfddm ! vertical physics: double diffusion 27 USE phycst ! physical constants 28 USE in_out_manager ! I/O manager 29 USE iom ! I/O manager library 30 USE lib_mpp ! MPP library 31 USE wrk_nemo ! Memory allocation 22 32 23 33 IMPLICIT NONE 24 34 PRIVATE 25 35 26 PUBLIC trd_tra ! called by all traXX modules 27 28 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: trdtx, trdty, trdt !: 36 PUBLIC trd_tra ! called by all tra_... modules 37 38 REAL(wp) :: r2dt ! time-step, = 2 rdttra except at nit000 (=rdttra) if neuler=0 39 40 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: trdtx, trdty, trdt ! use to store the temperature trends 29 41 30 42 !! * Substitutions 31 43 # include "domzgr_substitute.h90" 44 # include "zdfddm_substitute.h90" 32 45 # include "vectopt_loop_substitute.h90" 33 46 !!---------------------------------------------------------------------- 34 !! NEMO/OPA 4.0 , NEMO Consortium (2011)47 !! NEMO/OPA 3.3 , NEMO Consortium (2010) 35 48 !! $Id$ 36 49 !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) … … 39 52 40 53 INTEGER FUNCTION trd_tra_alloc() 41 !!--------------------------------------------------------------------- -------54 !!--------------------------------------------------------------------- 42 55 !! *** FUNCTION trd_tra_alloc *** 43 !!--------------------------------------------------------------------- -------56 !!--------------------------------------------------------------------- 44 57 ALLOCATE( trdtx(jpi,jpj,jpk) , trdty(jpi,jpj,jpk) , trdt(jpi,jpj,jpk) , STAT= trd_tra_alloc ) 45 58 ! … … 53 66 !! *** ROUTINE trd_tra *** 54 67 !! 55 !! ** Purpose : Dispatch all trends computation, e.g. vorticity, mld or 56 !! integral constraints 68 !! ** Purpose : pre-process tracer trends 57 69 !! 58 !! ** Method /usage : For the mixed-layer trend, the control surface can be either59 !! a mixed layer depth (time varying) or a fixed surface (jk level or bowl).60 !! Choose control surface with nn_ctls in namelist NAMTRD :61 !! nn_ctls = 0 : use mixed layer with density criterion62 !! nn_ctls = 1 : read index from file 'ctlsurf_idx'63 !! nn_ctls > 1 : use fixed level surface jk = nn_ctls64 !!---------------------------------------------------------------------- 65 !66 INTEGER , INTENT(in) :: kt ! time step67 CHARACTER(len=3) , INTENT(in) :: ctype ! tracers trends type 'TRA'/'TRC'68 INTEGER , INTENT(in) :: ktra ! tracerindex69 INTEGER , INTENT(in) :: ktrd ! tracer trend index70 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: ptrd ! tracer trend or flux71 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: pun ! velocity72 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: ptra ! Tracer variablea73 !74 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztrds75 !!---------------------------------------------------------------------- 76 70 !! ** Method : - mask the trend 71 !! - advection (ptra present) converte the incoming flux (U.T) 72 !! into trend (U.T => -U.grat(T)=div(U.T)-T.div(U)) through a 73 !! call to trd_tra_adv 74 !! - 'TRA' case : regroup T & S trends 75 !! - send the trends to trd_tra_mng (trdtrc) for further processing 76 !!---------------------------------------------------------------------- 77 INTEGER , INTENT(in) :: kt ! time step 78 CHARACTER(len=3) , INTENT(in) :: ctype ! tracers trends type 'TRA'/'TRC' 79 INTEGER , INTENT(in) :: ktra ! tracer index 80 INTEGER , INTENT(in) :: ktrd ! tracer trend index 81 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: ptrd ! tracer trend or flux 82 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: pun ! now velocity 83 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: ptra ! now tracer variable 84 ! 85 INTEGER :: jk ! loop indices 86 REAL(wp), POINTER, DIMENSION(:,:,:) :: zwt, zws, ztrdt, ztrds ! 3D workspace 87 !!---------------------------------------------------------------------- 88 ! 77 89 CALL wrk_alloc( jpi, jpj, jpk, ztrds ) 78 79 IF( .NOT. ALLOCATED( trdtx ) ) THEN 90 ! 91 IF( .NOT. ALLOCATED( trdtx ) ) THEN ! allocate trdtra arrays 80 92 IF( trd_tra_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'trd_tra : unable to allocate arrays' ) 81 93 ENDIF 82 83 ! Control of optional arguments 84 IF( ctype == 'TRA' .AND. ktra == jp_tem ) THEN 85 IF( PRESENT( ptra ) ) THEN 86 SELECT CASE( ktrd ) ! shift depending on the direction 87 CASE( jptra_trd_xad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'X', trdtx ) 88 CASE( jptra_trd_yad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'Y', trdty ) 89 CASE( jptra_trd_zad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'Z', trdt ) 90 END SELECT 91 ELSE 92 trdt(:,:,:) = ptrd(:,:,:) 93 IF( ktrd == jptra_trd_bbc .OR. ktrd == jptra_trd_qsr ) THEN 94 ztrds(:,:,:) = 0. 95 CALL trd_mod( trdt, ztrds, ktrd, ctype, kt ) 96 END IF 97 END IF 98 END IF 99 100 IF( ctype == 'TRA' .AND. ktra == jp_sal ) THEN 101 IF( PRESENT( ptra ) ) THEN 102 SELECT CASE( ktrd ) ! shift depending on the direction 103 CASE( jptra_trd_xad ) 104 CALL trd_tra_adv( ptrd, pun, ptra, 'X', ztrds ) 105 CALL trd_mod( trdtx, ztrds, ktrd, ctype, kt ) 106 CASE( jptra_trd_yad ) 107 CALL trd_tra_adv( ptrd, pun, ptra, 'Y', ztrds ) 108 CALL trd_mod( trdty, ztrds, ktrd, ctype, kt ) 109 CASE( jptra_trd_zad ) 110 CALL trd_tra_adv( ptrd, pun, ptra, 'Z', ztrds ) 111 CALL trd_mod( trdt , ztrds, ktrd, ctype, kt ) 112 END SELECT 113 ELSE 114 ztrds(:,:,:) = ptrd(:,:,:) 115 CALL trd_mod( trdt, ztrds, ktrd, ctype, kt ) 116 END IF 117 END IF 118 119 IF( ctype == 'TRC' ) THEN 120 ! 121 IF( PRESENT( ptra ) ) THEN 122 SELECT CASE( ktrd ) ! shift depending on the direction 123 CASE( jptra_trd_xad ) 124 CALL trd_tra_adv( ptrd, pun, ptra, 'X', ztrds ) 125 CALL trd_mod_trc( ztrds, ktra, ktrd, kt ) 126 CASE( jptra_trd_yad ) 127 CALL trd_tra_adv( ptrd, pun, ptra, 'Y', ztrds ) 128 CALL trd_mod_trc( ztrds, ktra, ktrd, kt ) 129 CASE( jptra_trd_zad ) 130 CALL trd_tra_adv( ptrd, pun, ptra, 'Z', ztrds ) 131 CALL trd_mod_trc( ztrds, ktra, ktrd, kt ) 132 END SELECT 133 ELSE 134 ztrds(:,:,:) = ptrd(:,:,:) 135 CALL trd_mod_trc( ztrds, ktra, ktrd, kt ) 136 END IF 94 95 IF( ctype == 'TRA' .AND. ktra == jp_tem ) THEN !== Temperature trend ==! 96 ! 97 SELECT CASE( ktrd ) 98 ! ! advection: transform the advective flux into a trend 99 CASE( jptra_xad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'X', trdtx ) 100 CASE( jptra_yad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'Y', trdty ) 101 CASE( jptra_zad ) ; CALL trd_tra_adv( ptrd, pun, ptra, 'Z', trdt ) 102 CASE( jptra_bbc, & ! qsr, bbc: on temperature only, send to trd_tra_mng 103 & jptra_qsr ) ; trdt(:,:,:) = ptrd(:,:,:) * tmask(:,:,:) 104 ztrds(:,:,:) = 0._wp 105 CALL trd_tra_mng( trdt, ztrds, ktrd, kt ) 106 CASE DEFAULT ! other trends: masked trends 107 trdt(:,:,:) = ptrd(:,:,:) * tmask(:,:,:) ! mask & store 108 END SELECT 109 ! 110 ENDIF 111 112 IF( ctype == 'TRA' .AND. ktra == jp_sal ) THEN !== Salinity trends ==! 113 ! 114 SELECT CASE( ktrd ) 115 ! ! advection: transform the advective flux into a trend 116 ! ! and send T & S trends to trd_tra_mng 117 CASE( jptra_xad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'X' , ztrds ) 118 CALL trd_tra_mng( trdtx, ztrds, ktrd, kt ) 119 CASE( jptra_yad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'Y' , ztrds ) 120 CALL trd_tra_mng( trdty, ztrds, ktrd, kt ) 121 CASE( jptra_zad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'Z' , ztrds ) 122 CALL trd_tra_mng( trdt , ztrds, ktrd, kt ) 123 CASE( jptra_zdfp ) ! diagnose the "PURE" Kz trend (here: just before the swap) 124 ! ! iso-neutral diffusion case otherwise jptra_zdf is "PURE" 125 CALL wrk_alloc( jpi, jpj, jpk, zwt, zws, ztrdt ) 126 ! 127 zwt(:,:, 1 ) = 0._wp ; zws(:,:, 1 ) = 0._wp ! vertical diffusive fluxes 128 zwt(:,:,jpk) = 0._wp ; zws(:,:,jpk) = 0._wp 129 DO jk = 2, jpk 130 zwt(:,:,jk) = avt(:,:,jk) * ( tsa(:,:,jk-1,jp_tem) - tsa(:,:,jk,jp_tem) ) / fse3w(:,:,jk) * tmask(:,:,jk) 131 zws(:,:,jk) = fsavs(:,:,jk) * ( tsa(:,:,jk-1,jp_sal) - tsa(:,:,jk,jp_sal) ) / fse3w(:,:,jk) * tmask(:,:,jk) 132 END DO 133 ! 134 ztrdt(:,:,jpk) = 0._wp ; ztrds(:,:,jpk) = 0._wp 135 DO jk = 1, jpkm1 136 ztrdt(:,:,jk) = ( zwt(:,:,jk) - zwt(:,:,jk+1) ) / fse3t(:,:,jk) 137 ztrds(:,:,jk) = ( zws(:,:,jk) - zws(:,:,jk+1) ) / fse3t(:,:,jk) 138 END DO 139 CALL trd_tra_mng( ztrdt, ztrds, jptra_zdfp, kt ) 140 ! 141 CALL wrk_dealloc( jpi, jpj, jpk, zwt, zws, ztrdt ) 142 ! 143 CASE DEFAULT ! other trends: mask and send T & S trends to trd_tra_mng 144 ztrds(:,:,:) = ptrd(:,:,:) * tmask(:,:,:) 145 CALL trd_tra_mng( trdt, ztrds, ktrd, kt ) 146 END SELECT 147 ENDIF 148 149 IF( ctype == 'TRC' ) THEN !== passive tracer trend ==! 150 ! 151 SELECT CASE( ktrd ) 152 ! ! advection: transform the advective flux into a masked trend 153 CASE( jptra_xad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'X', ztrds ) 154 CASE( jptra_yad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'Y', ztrds ) 155 CASE( jptra_zad ) ; CALL trd_tra_adv( ptrd , pun , ptra, 'Z', ztrds ) 156 CASE DEFAULT ! other trends: just masked 157 ztrds(:,:,:) = ptrd(:,:,:) * tmask(:,:,:) 158 END SELECT 159 ! ! send trend to trd_trc 160 CALL trd_trc( ztrds, ktra, ktrd, kt ) 137 161 ! 138 162 ENDIF … … 147 171 !! *** ROUTINE trd_tra_adv *** 148 172 !! 149 !! ** Purpose : transformed the i-, j- or k-advective flux into thes 150 !! i-, j- or k-advective trends, resp. 151 !! ** Method : i-advective trends = -un. di-1[T] = -( di-1[fi] - tn di-1[un] ) 152 !! k-advective trends = -un. di-1[T] = -( dj-1[fi] - tn dj-1[un] ) 153 !! k-advective trends = -un. di+1[T] = -( dk+1[fi] - tn dk+1[un] ) 154 !!---------------------------------------------------------------------- 155 REAL(wp) , INTENT(in ), DIMENSION(jpi,jpj,jpk) :: pf ! advective flux in one direction 156 REAL(wp) , INTENT(in ), DIMENSION(jpi,jpj,jpk) :: pun ! now velocity in one direction 157 REAL(wp) , INTENT(in ), DIMENSION(jpi,jpj,jpk) :: ptn ! now or before tracer 158 CHARACTER(len=1), INTENT(in ) :: cdir ! X/Y/Z direction 159 REAL(wp) , INTENT(out), DIMENSION(jpi,jpj,jpk) :: ptrd ! advective trend in one direction 173 !! ** Purpose : transformed a advective flux into a masked advective trends 174 !! 175 !! ** Method : use the following transformation: -div(U.T) = - U grad(T) + T.div(U) 176 !! i-advective trends = -un. di-1[T] = -( di-1[fi] - tn di-1[un] ) 177 !! j-advective trends = -un. di-1[T] = -( dj-1[fi] - tn dj-1[un] ) 178 !! k-advective trends = -un. di+1[T] = -( dk+1[fi] - tn dk+1[un] ) 179 !! where fi is the incoming advective flux. 180 !!---------------------------------------------------------------------- 181 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in ) :: pf ! advective flux in one direction 182 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in ) :: pun ! now velocity in one direction 183 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in ) :: ptn ! now or before tracer 184 CHARACTER(len=1) , INTENT(in ) :: cdir ! X/Y/Z direction 185 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( out) :: ptrd ! advective trend in one direction 160 186 ! 161 187 INTEGER :: ji, jj, jk ! dummy loop indices 162 INTEGER :: ii, ij, ik ! index shift function of the direction 163 REAL(wp) :: zbtr ! local scalar 164 !!---------------------------------------------------------------------- 165 166 SELECT CASE( cdir ) ! shift depending on the direction 167 CASE( 'X' ) ; ii = 1 ; ij = 0 ; ik = 0 ! i-advective trend 168 CASE( 'Y' ) ; ii = 0 ; ij = 1 ; ik = 0 ! j-advective trend 169 CASE( 'Z' ) ; ii = 0 ; ij = 0 ; ik =-1 ! k-advective trend 188 INTEGER :: ii, ij, ik ! index shift as function of the direction 189 !!---------------------------------------------------------------------- 190 ! 191 SELECT CASE( cdir ) ! shift depending on the direction 192 CASE( 'X' ) ; ii = 1 ; ij = 0 ; ik = 0 ! i-trend 193 CASE( 'Y' ) ; ii = 0 ; ij = 1 ; ik = 0 ! j-trend 194 CASE( 'Z' ) ; ii = 0 ; ij = 0 ; ik =-1 ! k-trend 170 195 END SELECT 171 172 ! ! set to zero uncomputed values 173 ptrd(jpi,:,:) = 0.e0 ; ptrd(1,:,:) = 0.e0 174 ptrd(:,jpj,:) = 0.e0 ; ptrd(:,1,:) = 0.e0 175 ptrd(:,:,jpk) = 0.e0 176 ! 177 ! 178 DO jk = 1, jpkm1 196 ! 197 ! ! set to zero uncomputed values 198 ptrd(jpi,:,:) = 0._wp ; ptrd(1,:,:) = 0._wp 199 ptrd(:,jpj,:) = 0._wp ; ptrd(:,1,:) = 0._wp 200 ptrd(:,:,jpk) = 0._wp 201 ! 202 DO jk = 1, jpkm1 ! advective trend 179 203 DO jj = 2, jpjm1 180 204 DO ji = fs_2, fs_jpim1 ! vector opt. 181 zbtr = 1.e0/ ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )182 ptrd(ji,jj,jk) = - zbtr * ( pf (ji,jj,jk) - pf (ji-ii,jj-ij,jk-ik)&183 & - ( pun(ji,jj,jk) - pun(ji-ii,jj-ij,jk-ik) ) * ptn(ji,jj,jk))205 ptrd(ji,jj,jk) = - ( pf (ji,jj,jk) - pf (ji-ii,jj-ij,jk-ik) & 206 & - ( pun(ji,jj,jk) - pun(ji-ii,jj-ij,jk-ik) ) * ptn(ji,jj,jk) ) & 207 & / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) * tmask(ji,jj,jk) 184 208 END DO 185 209 END DO … … 188 212 END SUBROUTINE trd_tra_adv 189 213 190 # else 191 !!---------------------------------------------------------------------- 192 !! Default case : Dummy module No trend diagnostics 193 !!---------------------------------------------------------------------- 194 USE par_oce ! ocean variables trends 195 CONTAINS 196 SUBROUTINE trd_tra( kt, ctype, ktra, ktrd, ptrd, pu, ptra ) 197 !!---------------------------------------------------------------------- 198 INTEGER , INTENT(in) :: kt ! time step 199 CHARACTER(len=3) , INTENT(in) :: ctype ! tracers trends type 'TRA'/'TRC' 200 INTEGER , INTENT(in) :: ktra ! tracer index 201 INTEGER , INTENT(in) :: ktrd ! tracer trend index 202 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: ptrd ! tracer trend 203 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: pu ! velocity 204 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in), OPTIONAL :: ptra ! Tracer variable 205 WRITE(*,*) 'trd_3d: You should not have seen this print! error ?', ptrd(1,1,1), ptra(1,1,1), pu(1,1,1), & 206 & ktrd, ktra, ctype, kt 207 END SUBROUTINE trd_tra 208 # endif 214 215 SUBROUTINE trd_tra_mng( ptrdx, ptrdy, ktrd, kt ) 216 !!--------------------------------------------------------------------- 217 !! *** ROUTINE trd_tra_mng *** 218 !! 219 !! ** Purpose : Dispatch all tracer trends computation, e.g. 3D output, 220 !! integral constraints, potential energy, and/or 221 !! mixed layer budget. 222 !!---------------------------------------------------------------------- 223 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdx ! Temperature or U trend 224 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdy ! Salinity or V trend 225 INTEGER , INTENT(in ) :: ktrd ! tracer trend index 226 INTEGER , INTENT(in ) :: kt ! time step 227 !!---------------------------------------------------------------------- 228 229 IF( neuler == 0 .AND. kt == nit000 ) THEN ; r2dt = rdt ! = rdtra (restart with Euler time stepping) 230 ELSEIF( kt <= nit000 + 1) THEN ; r2dt = 2. * rdt ! = 2 rdttra (leapfrog) 231 ENDIF 232 233 ! ! 3D output of tracers trends using IOM interface 234 IF( ln_tra_trd ) CALL trd_tra_iom ( ptrdx, ptrdy, ktrd, kt ) 235 236 ! ! Integral Constraints Properties for tracers trends !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 237 IF( ln_glo_trd ) CALL trd_glo( ptrdx, ptrdy, ktrd, 'TRA', kt ) 238 239 ! ! Potential ENergy trends 240 IF( ln_PE_trd ) CALL trd_pen( ptrdx, ptrdy, ktrd, kt, r2dt ) 241 242 ! ! Mixed layer trends for active tracers 243 IF( ln_tra_mxl ) THEN 244 !----------------------------------------------------------------------------------------------- 245 ! W.A.R.N.I.N.G : 246 ! jptra_ldf : called by traldf.F90 247 ! at this stage we store: 248 ! - the lateral geopotential diffusion (here, lateral = horizontal) 249 ! - and the iso-neutral diffusion if activated 250 ! jptra_zdf : called by trazdf.F90 251 ! * in case of iso-neutral diffusion we store the vertical diffusion component in the 252 ! lateral trend including the K_z contrib, which will be removed later (see trd_mxl) 253 !----------------------------------------------------------------------------------------------- 254 255 SELECT CASE ( ktrd ) 256 CASE ( jptra_xad ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_xad, '3D' ) ! zonal advection 257 CASE ( jptra_yad ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_yad, '3D' ) ! merid. advection 258 CASE ( jptra_zad ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_zad, '3D' ) ! vertical advection 259 CASE ( jptra_ldf ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_ldf, '3D' ) ! lateral diffusion 260 CASE ( jptra_bbl ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_bbl, '3D' ) ! bottom boundary layer 261 CASE ( jptra_zdf ) 262 IF( ln_traldf_iso ) THEN ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_ldf, '3D' ) ! lateral diffusion (K_z) 263 ELSE ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_zdf, '3D' ) ! vertical diffusion (K_z) 264 ENDIF 265 CASE ( jptra_dmp ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_dmp, '3D' ) ! internal 3D restoring (tradmp) 266 CASE ( jptra_qsr ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_for, '3D' ) ! air-sea : penetrative sol radiat 267 CASE ( jptra_nsr ) ; ptrdx(:,:,2:jpk) = 0._wp ; ptrdy(:,:,2:jpk) = 0._wp 268 CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_for, '2D' ) ! air-sea : non penetr sol radiation 269 CASE ( jptra_bbc ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_bbc, '3D' ) ! bottom bound cond (geoth flux) 270 CASE ( jptra_npc ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_npc, '3D' ) ! non penetr convect adjustment 271 CASE ( jptra_atf ) ; CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_atf, '3D' ) ! asselin time filter (last trend) 272 ! 273 CALL trd_mxl( kt, r2dt ) ! trends: Mixed-layer (output) 274 END SELECT 275 ! 276 ENDIF 277 ! 278 END SUBROUTINE trd_tra_mng 279 280 281 SUBROUTINE trd_tra_iom( ptrdx, ptrdy, ktrd, kt ) 282 !!--------------------------------------------------------------------- 283 !! *** ROUTINE trd_tra_iom *** 284 !! 285 !! ** Purpose : output 3D tracer trends using IOM 286 !!---------------------------------------------------------------------- 287 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdx ! Temperature or U trend 288 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdy ! Salinity or V trend 289 INTEGER , INTENT(in ) :: ktrd ! tracer trend index 290 INTEGER , INTENT(in ) :: kt ! time step 291 !! 292 INTEGER :: ji, jj, jk ! dummy loop indices 293 INTEGER :: ikbu, ikbv ! local integers 294 REAL(wp), POINTER, DIMENSION(:,:) :: z2dx, z2dy ! 2D workspace 295 !!---------------------------------------------------------------------- 296 ! 297 !!gm Rq: mask the trends already masked in trd_tra, but lbc_lnk should probably be added 298 ! 299 SELECT CASE( ktrd ) 300 CASE( jptra_xad ) ; CALL iom_put( "ttrd_xad" , ptrdx ) ! x- horizontal advection 301 CALL iom_put( "strd_xad" , ptrdy ) 302 CASE( jptra_yad ) ; CALL iom_put( "ttrd_yad" , ptrdx ) ! y- horizontal advection 303 CALL iom_put( "strd_yad" , ptrdy ) 304 CASE( jptra_zad ) ; CALL iom_put( "ttrd_zad" , ptrdx ) ! z- vertical advection 305 CALL iom_put( "strd_zad" , ptrdy ) 306 IF( .NOT. lk_vvl ) THEN ! cst volume : adv flux through z=0 surface 307 CALL wrk_alloc( jpi, jpj, z2dx, z2dy ) 308 z2dx(:,:) = wn(:,:,1) * tsn(:,:,1,jp_tem) / fse3t(:,:,1) 309 z2dy(:,:) = wn(:,:,1) * tsn(:,:,1,jp_sal) / fse3t(:,:,1) 310 CALL iom_put( "ttrd_sad", z2dx ) 311 CALL iom_put( "strd_sad", z2dy ) 312 CALL wrk_dealloc( jpi, jpj, z2dx, z2dy ) 313 ENDIF 314 CASE( jptra_ldf ) ; CALL iom_put( "ttrd_ldf" , ptrdx ) ! lateral diffusion 315 CALL iom_put( "strd_ldf" , ptrdy ) 316 CASE( jptra_zdf ) ; CALL iom_put( "ttrd_zdf" , ptrdx ) ! vertical diffusion (including Kz contribution) 317 CALL iom_put( "strd_zdf" , ptrdy ) 318 CASE( jptra_zdfp ) ; CALL iom_put( "ttrd_zdfp", ptrdx ) ! PURE vertical diffusion (no isoneutral contribution) 319 CALL iom_put( "strd_zdfp", ptrdy ) 320 CASE( jptra_dmp ) ; CALL iom_put( "ttrd_dmp" , ptrdx ) ! internal restoring (damping) 321 CALL iom_put( "strd_dmp" , ptrdy ) 322 CASE( jptra_bbl ) ; CALL iom_put( "ttrd_bbl" , ptrdx ) ! bottom boundary layer 323 CALL iom_put( "strd_bbl" , ptrdy ) 324 CASE( jptra_npc ) ; CALL iom_put( "ttrd_npc" , ptrdx ) ! static instability mixing 325 CALL iom_put( "strd_npc" , ptrdy ) 326 CASE( jptra_nsr ) ; CALL iom_put( "ttrd_qns" , ptrdx ) ! surface forcing + runoff (ln_rnf=T) 327 CALL iom_put( "strd_cdt" , ptrdy ) 328 CASE( jptra_qsr ) ; CALL iom_put( "ttrd_qsr" , ptrdx ) ! penetrative solar radiat. (only on temperature) 329 CASE( jptra_bbc ) ; CALL iom_put( "ttrd_bbc" , ptrdx ) ! geothermal heating (only on temperature) 330 CASE( jptra_atf ) ; CALL iom_put( "ttrd_atf" , ptrdx ) ! asselin time Filter 331 CALL iom_put( "strd_atf" , ptrdy ) 332 END SELECT 333 ! 334 END SUBROUTINE trd_tra_iom 335 209 336 !!====================================================================== 210 337 END MODULE trdtra -
trunk/NEMOGCM/NEMO/OPA_SRC/TRD/trdvor.F90
r3294 r4990 4 4 !! Ocean diagnostics: momentum trends 5 5 !!===================================================================== 6 !! History : 1.0 ! 04-2006 (L. Brunier, A-M. Treguier) Original code 7 !! 2.0 ! 04-2008 (C. Talandier) New trends organization 6 !! History : 1.0 ! 2006-01 (L. Brunier, A-M. Treguier) Original code 7 !! 2.0 ! 2008-04 (C. Talandier) New trends organization 8 !! 3.5 ! 2012-02 (G. Madec) regroup beta.V computation with pvo trend 8 9 !!---------------------------------------------------------------------- 9 #if defined key_trdvor || defined key_esopa 10 !!---------------------------------------------------------------------- 11 !! 'key_trdvor' : momentum trend diagnostics 10 12 11 !!---------------------------------------------------------------------- 13 12 !! trd_vor : momentum trends averaged over the depth … … 17 16 USE oce ! ocean dynamics and tracers variables 18 17 USE dom_oce ! ocean space and time domain variables 19 USE trd mod_oce ! ocean variables trends18 USE trd_oce ! trends: ocean variables 20 19 USE zdf_oce ! ocean vertical physics 21 USE in_out_manager ! I/O manager20 USE sbc_oce ! surface boundary condition: ocean 22 21 USE phycst ! Define parameters for the routines 23 22 USE ldfdyn_oce ! ocean active tracers: lateral physics 24 23 USE dianam ! build the name of file (routine) 25 24 USE zdfmxl ! mixed layer depth 25 USE lbclnk ! ocean lateral boundary conditions (or mpp link) 26 USE in_out_manager ! I/O manager 26 27 USE ioipsl ! NetCDF library 27 USE lbclnk ! ocean lateral boundary conditions (or mpp link)28 28 USE lib_mpp ! MPP library 29 29 USE wrk_nemo ! Memory allocation 30 31 30 32 31 IMPLICIT NONE … … 37 36 END INTERFACE 38 37 39 PUBLIC trd_vor ! routine called by step.F90 40 PUBLIC trd_vor_zint ! routine called by dynamics routines 38 PUBLIC trd_vor ! routine called by trddyn.F90 41 39 PUBLIC trd_vor_init ! routine called by opa.F90 42 40 PUBLIC trd_vor_alloc ! routine called by nemogcm.F90 … … 80 78 IF( trd_vor_alloc /= 0 ) CALL ctl_warn('trd_vor_alloc: failed to allocate arrays') 81 79 END FUNCTION trd_vor_alloc 80 81 82 SUBROUTINE trd_vor( putrd, pvtrd, ktrd, kt ) 83 !!---------------------------------------------------------------------- 84 !! *** ROUTINE trd_vor *** 85 !! 86 !! ** Purpose : computation of cumulated trends over analysis period 87 !! and make outputs (NetCDF or DIMG format) 88 !!---------------------------------------------------------------------- 89 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: putrd, pvtrd ! U and V trends 90 INTEGER , INTENT(in ) :: ktrd ! trend index 91 INTEGER , INTENT(in ) :: kt ! time step 92 ! 93 INTEGER :: ji, jj ! dummy loop indices 94 REAL(wp), POINTER, DIMENSION(:,:) :: ztswu, ztswv ! 2D workspace 95 !!---------------------------------------------------------------------- 96 97 CALL wrk_alloc( jpi, jpj, ztswu, ztswv ) 98 99 SELECT CASE( ktrd ) 100 CASE( jpdyn_hpg ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_prg ) ! Hydrostatique Pressure Gradient 101 CASE( jpdyn_keg ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_keg ) ! KE Gradient 102 CASE( jpdyn_rvo ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_rvo ) ! Relative Vorticity 103 CASE( jpdyn_pvo ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_pvo ) ! Planetary Vorticity Term 104 CASE( jpdyn_ldf ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_ldf ) ! Horizontal Diffusion 105 CASE( jpdyn_zad ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_zad ) ! Vertical Advection 106 CASE( jpdyn_spg ) ; CALL trd_vor_zint( putrd, pvtrd, jpvor_spg ) ! Surface Pressure Grad. 107 CASE( jpdyn_zdf ) ! Vertical Diffusion 108 ztswu(:,:) = 0.e0 ; ztswv(:,:) = 0.e0 109 DO jj = 2, jpjm1 ! wind stress trends 110 DO ji = fs_2, fs_jpim1 ! vector opt. 111 ztswu(ji,jj) = 0.5 * ( utau_b(ji,jj) + utau(ji,jj) ) / ( fse3u(ji,jj,1) * rau0 ) 112 ztswv(ji,jj) = 0.5 * ( vtau_b(ji,jj) + vtau(ji,jj) ) / ( fse3v(ji,jj,1) * rau0 ) 113 END DO 114 END DO 115 ! 116 CALL trd_vor_zint( putrd, pvtrd, jpvor_zdf ) ! zdf trend including surf./bot. stresses 117 CALL trd_vor_zint( ztswu, ztswv, jpvor_swf ) ! surface wind stress 118 CASE( jpdyn_bfr ) 119 CALL trd_vor_zint( putrd, pvtrd, jpvor_bfr ) ! Bottom stress 120 ! 121 CASE( jpdyn_atf ) ! last trends: perform the output of 2D vorticity trends 122 CALL trd_vor_iom( kt ) 123 END SELECT 124 ! 125 CALL wrk_dealloc( jpi, jpj, ztswu, ztswv ) 126 ! 127 END SUBROUTINE trd_vor 82 128 83 129 … … 109 155 !! trends output in netCDF format using ioipsl 110 156 !!---------------------------------------------------------------------- 111 !112 157 INTEGER , INTENT(in ) :: ktrd ! ocean trend index 113 158 REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: putrdvor ! u vorticity trend … … 131 176 ! ===================================== 132 177 133 SELECT CASE (ktrd)134 ! 135 CASE (jpvor_bfr) ! bottom friction178 SELECT CASE( ktrd ) 179 ! 180 CASE( jpvor_bfr ) ! bottom friction 136 181 DO jj = 2, jpjm1 137 182 DO ji = fs_2, fs_jpim1 … … 143 188 END DO 144 189 ! 145 CASE (jpvor_swf) ! wind stress190 CASE( jpvor_swf ) ! wind stress 146 191 zudpvor(:,:) = putrdvor(:,:) * fse3u(:,:,1) * e1u(:,:) * umask(:,:,1) 147 192 zvdpvor(:,:) = pvtrdvor(:,:) * fse3v(:,:,1) * e2v(:,:) * vmask(:,:,1) … … 154 199 155 200 ! Curl 156 DO ji =1,jpim1157 DO jj =1,jpjm1201 DO ji = 1, jpim1 202 DO jj = 1, jpjm1 158 203 vortrd(ji,jj,ktrd) = ( zvdpvor(ji+1,jj) - zvdpvor(ji,jj) & 159 204 & - ( zudpvor(ji,jj+1) - zudpvor(ji,jj) ) ) / ( e1f(ji,jj) * e2f(ji,jj) ) … … 229 274 END DO 230 275 231 ! Save Beta.V term to avoid average before Curl232 ! Beta.V : intergration, noaverage233 IF( ktrd == jpvor_ bev) THEN276 ! Planetary vorticity: 2nd computation (Beta.V term) store the vertical sum 277 ! as Beta.V term need intergration, not average 278 IF( ktrd == jpvor_pvo ) THEN 234 279 zubet(:,:) = zudpvor(:,:) 235 280 zvbet(:,:) = zvdpvor(:,:) 236 ENDIF 237 238 ! Average except for Beta.V 281 DO ji = 1, jpim1 282 DO jj = 1, jpjm1 283 vortrd(ji,jj,jpvor_bev) = ( zvbet(ji+1,jj) - zvbet(ji,jj) & 284 & - ( zubet(ji,jj+1) - zubet(ji,jj) ) ) / ( e1f(ji,jj) * e2f(ji,jj) ) 285 END DO 286 END DO 287 ! Average of the Curl and Surface mask 288 vortrd(:,:,jpvor_bev) = vortrd(:,:,jpvor_bev) * hur(:,:) * fmask(:,:,1) 289 ENDIF 290 ! 291 ! Average 239 292 zudpvor(:,:) = zudpvor(:,:) * hur(:,:) 240 293 zvdpvor(:,:) = zvdpvor(:,:) * hvr(:,:) 241 294 ! 242 295 ! Curl 243 296 DO ji=1,jpim1 … … 247 300 END DO 248 301 END DO 249 250 302 ! Surface mask 251 303 vortrd(:,:,ktrd) = vortrd(:,:,ktrd) * fmask(:,:,1) 252 253 ! Special treatement for the Beta.V term254 ! Compute the Curl of the Beta.V term which is not averaged255 IF( ktrd == jpvor_bev ) THEN256 DO ji=1,jpim1257 DO jj=1,jpjm1258 vortrd(ji,jj,jpvor_bev) = ( zvbet(ji+1,jj) - zvbet(ji,jj) &259 & - ( zubet(ji,jj+1) - zubet(ji,jj) ) ) / ( e1f(ji,jj) * e2f(ji,jj) )260 END DO261 END DO262 263 ! Average on the Curl264 vortrd(:,:,jpvor_bev) = vortrd(:,:,jpvor_bev) * hur(:,:)265 266 ! Surface mask267 vortrd(:,:,jpvor_bev) = vortrd(:,:,jpvor_bev) * fmask(:,:,1)268 ENDIF269 304 270 305 IF( ndebug /= 0 ) THEN … … 278 313 279 314 280 SUBROUTINE trd_vor ( kt )315 SUBROUTINE trd_vor_iom( kt ) 281 316 !!---------------------------------------------------------------------- 282 317 !! *** ROUTINE trd_vor *** … … 285 320 !! and make outputs (NetCDF or DIMG format) 286 321 !!---------------------------------------------------------------------- 287 ! 288 INTEGER, INTENT(in) :: kt ! ocean time-step index 322 INTEGER , INTENT(in ) :: kt ! time step 289 323 ! 290 324 INTEGER :: ji, jj, jk, jl ! dummy loop indices … … 305 339 306 340 IF( kt > nit000 ) vor_avrb(:,:) = vor_avr(:,:) 307 308 IF( ndebug /= 0 ) THEN309 WRITE(numout,*) ' debuging trd_vor: I.1 done '310 CALL FLUSH(numout)311 ENDIF312 341 313 342 ! I.2 vertically integrated vorticity … … 330 359 331 360 ! Curl 332 DO ji =1,jpim1333 DO jj =1,jpjm1361 DO ji = 1, jpim1 362 DO jj = 1, jpjm1 334 363 vor_avr(ji,jj) = ( ( zvn(ji+1,jj) - zvn(ji,jj) ) & 335 364 & - ( zun(ji,jj+1) - zun(ji,jj) ) ) / ( e1f(ji,jj) * e2f(ji,jj) ) * fmask(ji,jj,1) … … 337 366 END DO 338 367 339 IF( ndebug /= 0 ) THEN340 WRITE(numout,*) ' debuging trd_vor: I.2 done'341 CALL FLUSH(numout)342 ENDIF343 344 368 ! ================================= 345 369 ! II. Cumulated trends … … 351 375 vor_avrbb(:,:) = vor_avrb(:,:) 352 376 vor_avrbn(:,:) = vor_avr (:,:) 353 ENDIF354 355 IF( ndebug /= 0 ) THEN356 WRITE(numout,*) ' debuging trd_vor: I1.1 done'357 CALL FLUSH(numout)358 377 ENDIF 359 378 … … 371 390 ENDIF 372 391 373 IF( ndebug /= 0 ) THEN374 WRITE(numout,*) ' debuging trd_vor: II.2 done'375 CALL FLUSH(numout)376 ENDIF377 378 392 ! ============================================= 379 393 ! III. Output in netCDF + residual computation … … 391 405 vor_avrtot(:,:) = ( vor_avr(:,:) - vor_avrbn(:,:) + vor_avrb(:,:) - vor_avrbb(:,:) ) * zmean 392 406 393 IF( ndebug /= 0 ) THEN394 WRITE(numout,*) ' zmean = ',zmean395 WRITE(numout,*) ' debuging trd_vor: III.1 done'396 CALL FLUSH(numout)397 ENDIF398 407 399 408 ! III.2 compute residual … … 406 415 CALL lbc_lnk( vor_avrres, 'F', 1. ) 407 416 408 IF( ndebug /= 0 ) THEN409 WRITE(numout,*) ' debuging trd_vor: III.2 done'410 CALL FLUSH(numout)411 ENDIF412 417 413 418 ! III.3 time evolution array swap … … 415 420 vor_avrbb(:,:) = vor_avrb(:,:) 416 421 vor_avrbn(:,:) = vor_avr (:,:) 417 418 IF( ndebug /= 0 ) THEN419 WRITE(numout,*) ' debuging trd_vor: III.3 done'420 CALL FLUSH(numout)421 ENDIF422 422 ! 423 423 nmoydpvor = 0 … … 463 463 CALL wrk_dealloc( jpi, jpj, zun, zvn ) 464 464 ! 465 END SUBROUTINE trd_vor 465 END SUBROUTINE trd_vor_iom 466 466 467 467 … … 587 587 END SUBROUTINE trd_vor_init 588 588 589 #else590 !!----------------------------------------------------------------------591 !! Default option : Empty module592 !!----------------------------------------------------------------------593 INTERFACE trd_vor_zint594 MODULE PROCEDURE trd_vor_zint_2d, trd_vor_zint_3d595 END INTERFACE596 CONTAINS597 SUBROUTINE trd_vor( kt ) ! Empty routine598 WRITE(*,*) 'trd_vor: You should not have seen this print! error?', kt599 END SUBROUTINE trd_vor600 SUBROUTINE trd_vor_zint_2d( putrdvor, pvtrdvor, ktrd )601 REAL, DIMENSION(:,:), INTENT( inout ) :: putrdvor, pvtrdvor602 INTEGER, INTENT( in ) :: ktrd ! ocean trend index603 WRITE(*,*) 'trd_vor_zint_2d: You should not have seen this print! error?', putrdvor(1,1), pvtrdvor(1,1), ktrd604 END SUBROUTINE trd_vor_zint_2d605 SUBROUTINE trd_vor_zint_3d( putrdvor, pvtrdvor, ktrd )606 REAL, DIMENSION(:,:,:), INTENT( inout ) :: putrdvor, pvtrdvor607 INTEGER, INTENT( in ) :: ktrd ! ocean trend index608 WRITE(*,*) 'trd_vor_zint_3d: You should not have seen this print! error?', putrdvor(1,1,1), pvtrdvor(1,1,1), ktrd609 END SUBROUTINE trd_vor_zint_3d610 SUBROUTINE trd_vor_init ! Empty routine611 WRITE(*,*) 'trd_vor_init: You should not have seen this print! error?'612 END SUBROUTINE trd_vor_init613 #endif614 589 !!====================================================================== 615 590 END MODULE trdvor -
trunk/NEMOGCM/NEMO/OPA_SRC/TRD/trdvor_oce.F90
r2715 r4990 4 4 !! Ocean trends : set vorticity trend variables 5 5 !!====================================================================== 6 !! History : 9.0 ! 04-2006 (L. Brunier, A-M. Treguier) Original code6 !! History : 1.0 ! 04-2006 (L. Brunier, A-M. Treguier) Original code 7 7 !!---------------------------------------------------------------------- 8 9 !!---------------------------------------------------------------------- 8 10 9 USE par_oce ! ocean parameters 11 10 … … 13 12 PRIVATE 14 13 15 #if defined key_trdvor16 LOGICAL, PUBLIC, PARAMETER :: lk_trdvor = .TRUE. !: momentum trend flag17 #else18 LOGICAL, PUBLIC, PARAMETER :: lk_trdvor = .FALSE. !: momentum trend flag19 #endif20 14 ! !!* vorticity trends index 21 15 INTEGER, PUBLIC, PARAMETER :: jpltot_vor = 11 !: Number of vorticity trend terms
Note: See TracChangeset
for help on using the changeset viewer.