[3] | 1 | MODULE divcur |
---|
| 2 | !!============================================================================== |
---|
| 3 | !! *** MODULE divcur *** |
---|
| 4 | !! Ocean diagnostic variable : horizontal divergence and relative vorticity |
---|
| 5 | !!============================================================================== |
---|
[2528] | 6 | !! History : OPA ! 1987-06 (P. Andrich, D. L Hostis) Original code |
---|
| 7 | !! 4.0 ! 1991-11 (G. Madec) |
---|
| 8 | !! 6.0 ! 1993-03 (M. Guyon) symetrical conditions |
---|
| 9 | !! 7.0 ! 1996-01 (G. Madec) s-coordinates |
---|
| 10 | !! 8.0 ! 1997-06 (G. Madec) lateral boundary cond., lbc |
---|
| 11 | !! 8.1 ! 1997-08 (J.M. Molines) Open boundaries |
---|
| 12 | !! 8.2 ! 2000-03 (G. Madec) no slip accurate |
---|
| 13 | !! NEMO 1.0 ! 2002-09 (G. Madec, E. Durand) Free form, F90 |
---|
| 14 | !! - ! 2005-01 (J. Chanut) Unstructured open boundaries |
---|
| 15 | !! - ! 2003-08 (G. Madec) merged of cur and div, free form, F90 |
---|
| 16 | !! - ! 2005-01 (J. Chanut, A. Sellar) unstructured open boundaries |
---|
| 17 | !! 3.3 ! 2010-09 (D.Storkey and E.O'Dea) bug fixes for BDY module |
---|
| 18 | !! - ! 2010-10 (R. Furner, G. Madec) runoff and cla added directly here |
---|
[5120] | 19 | !! 3.6 ! 2014-11 (P. Mathiot) isf added directly here |
---|
[2528] | 20 | !!---------------------------------------------------------------------- |
---|
[3] | 21 | |
---|
| 22 | !!---------------------------------------------------------------------- |
---|
| 23 | !! div_cur : Compute the horizontal divergence and relative |
---|
| 24 | !! vorticity fields |
---|
| 25 | !!---------------------------------------------------------------------- |
---|
| 26 | USE oce ! ocean dynamics and tracers |
---|
| 27 | USE dom_oce ! ocean space and time domain |
---|
[4990] | 28 | USE sbc_oce, ONLY : ln_rnf, nn_isf ! surface boundary condition: ocean |
---|
[2528] | 29 | USE sbcrnf ! river runoff |
---|
[4990] | 30 | USE sbcisf ! ice shelf |
---|
[2528] | 31 | USE cla ! cross land advection (cla_div routine) |
---|
[3] | 32 | USE in_out_manager ! I/O manager |
---|
| 33 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
[2715] | 34 | USE lib_mpp ! MPP library |
---|
[3294] | 35 | USE wrk_nemo ! Memory Allocation |
---|
| 36 | USE timing ! Timing |
---|
[7845] | 37 | USE iom ! I/O Manager for dyn_vrt_dia_2d |
---|
[3] | 38 | |
---|
| 39 | IMPLICIT NONE |
---|
| 40 | PRIVATE |
---|
| 41 | |
---|
[7845] | 42 | PUBLIC div_cur ! routine called by step.F90 and istate.F90 |
---|
| 43 | PUBLIC dyn_vrt_dia_3d ! routine called by various modules |
---|
| 44 | PUBLIC dyn_vrt_dia_2d ! routine called by various modules |
---|
[3] | 45 | |
---|
| 46 | !! * Substitutions |
---|
| 47 | # include "domzgr_substitute.h90" |
---|
| 48 | # include "vectopt_loop_substitute.h90" |
---|
| 49 | !!---------------------------------------------------------------------- |
---|
[2528] | 50 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
[6486] | 51 | !! $Id$ |
---|
[2528] | 52 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[3] | 53 | !!---------------------------------------------------------------------- |
---|
| 54 | CONTAINS |
---|
| 55 | |
---|
| 56 | #if defined key_noslip_accurate |
---|
| 57 | !!---------------------------------------------------------------------- |
---|
[2715] | 58 | !! 'key_noslip_accurate' 2nd order interior + 4th order at the coast |
---|
[3] | 59 | !!---------------------------------------------------------------------- |
---|
| 60 | |
---|
| 61 | SUBROUTINE div_cur( kt ) |
---|
| 62 | !!---------------------------------------------------------------------- |
---|
| 63 | !! *** ROUTINE div_cur *** |
---|
| 64 | !! |
---|
| 65 | !! ** Purpose : compute the horizontal divergence and the relative |
---|
[2715] | 66 | !! vorticity at before and now time-step |
---|
[3] | 67 | !! |
---|
[2528] | 68 | !! ** Method : I. divergence : |
---|
[3] | 69 | !! - save the divergence computed at the previous time-step |
---|
| 70 | !! (note that the Asselin filter has not been applied on hdivb) |
---|
| 71 | !! - compute the now divergence given by : |
---|
| 72 | !! hdivn = 1/(e1t*e2t*e3t) ( di[e2u*e3u un] + dj[e1v*e3v vn] ) |
---|
[4990] | 73 | !! correct hdiv with runoff inflow (div_rnf), ice shelf melting (div_isf) |
---|
| 74 | !! and cross land flow (div_cla) |
---|
[2528] | 75 | !! II. vorticity : |
---|
[3] | 76 | !! - save the curl computed at the previous time-step |
---|
| 77 | !! rotb = rotn |
---|
| 78 | !! (note that the Asselin time filter has not been applied to rotb) |
---|
| 79 | !! - compute the now curl in tensorial formalism: |
---|
| 80 | !! rotn = 1/(e1f*e2f) ( di[e2v vn] - dj[e1u un] ) |
---|
| 81 | !! - Coastal boundary condition: 'key_noslip_accurate' defined, |
---|
| 82 | !! the no-slip boundary condition is computed using Schchepetkin |
---|
| 83 | !! and O'Brien (1996) scheme (i.e. 4th order at the coast). |
---|
| 84 | !! For example, along east coast, the one-sided finite difference |
---|
| 85 | !! approximation used for di[v] is: |
---|
[2528] | 86 | !! di[e2v vn] = 1/(e1f*e2f) * ( (e2v vn)(i) + (e2v vn)(i-1) + (e2v vn)(i-2) ) |
---|
[3] | 87 | !! |
---|
| 88 | !! ** Action : - update hdivb, hdivn, the before & now hor. divergence |
---|
| 89 | !! - update rotb , rotn , the before & now rel. vorticity |
---|
| 90 | !!---------------------------------------------------------------------- |
---|
[2715] | 91 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
[2528] | 92 | ! |
---|
[2715] | 93 | INTEGER :: ji, jj, jk, jl ! dummy loop indices |
---|
| 94 | INTEGER :: ii, ij, ijt, iju, ierr ! local integer |
---|
| 95 | REAL(wp) :: zraur, zdep ! local scalar |
---|
[3294] | 96 | REAL(wp), POINTER, DIMENSION(:,:) :: zwu ! specific 2D workspace |
---|
| 97 | REAL(wp), POINTER, DIMENSION(:,:) :: zwv ! specific 2D workspace |
---|
[3] | 98 | !!---------------------------------------------------------------------- |
---|
[3294] | 99 | ! |
---|
| 100 | IF( nn_timing == 1 ) CALL timing_start('div_cur') |
---|
| 101 | ! |
---|
[6487] | 102 | CALL wrk_alloc( jpi , jpj+2, zwu ) |
---|
| 103 | CALL wrk_alloc( jpi+2, jpj , zwv ) |
---|
[3294] | 104 | ! |
---|
[3] | 105 | IF( kt == nit000 ) THEN |
---|
| 106 | IF(lwp) WRITE(numout,*) |
---|
| 107 | IF(lwp) WRITE(numout,*) 'div_cur : horizontal velocity divergence and relative vorticity' |
---|
| 108 | IF(lwp) WRITE(numout,*) '~~~~~~~ NOT optimal for auto-tasking case' |
---|
| 109 | ENDIF |
---|
| 110 | |
---|
| 111 | ! ! =============== |
---|
| 112 | DO jk = 1, jpkm1 ! Horizontal slab |
---|
| 113 | ! ! =============== |
---|
[2528] | 114 | ! |
---|
[3] | 115 | hdivb(:,:,jk) = hdivn(:,:,jk) ! time swap of div arrays |
---|
| 116 | rotb (:,:,jk) = rotn (:,:,jk) ! time swap of rot arrays |
---|
[2528] | 117 | ! |
---|
[3] | 118 | ! ! -------- |
---|
| 119 | ! Horizontal divergence ! div |
---|
| 120 | ! ! -------- |
---|
| 121 | DO jj = 2, jpjm1 |
---|
| 122 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 123 | hdivn(ji,jj,jk) = & |
---|
[455] | 124 | ( e2u(ji,jj)*fse3u(ji,jj,jk) * un(ji,jj,jk) - e2u(ji-1,jj )*fse3u(ji-1,jj ,jk) * un(ji-1,jj ,jk) & |
---|
| 125 | + e1v(ji,jj)*fse3v(ji,jj,jk) * vn(ji,jj,jk) - e1v(ji ,jj-1)*fse3v(ji ,jj-1,jk) * vn(ji ,jj-1,jk) ) & |
---|
[3] | 126 | / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) |
---|
| 127 | END DO |
---|
| 128 | END DO |
---|
| 129 | |
---|
[1792] | 130 | IF( .NOT. AGRIF_Root() ) THEN |
---|
| 131 | IF ((nbondi == 1).OR.(nbondi == 2)) hdivn(nlci-1 , : ,jk) = 0.e0 ! east |
---|
| 132 | IF ((nbondi == -1).OR.(nbondi == 2)) hdivn(2 , : ,jk) = 0.e0 ! west |
---|
| 133 | IF ((nbondj == 1).OR.(nbondj == 2)) hdivn(: ,nlcj-1 ,jk) = 0.e0 ! north |
---|
| 134 | IF ((nbondj == -1).OR.(nbondj == 2)) hdivn(: ,2 ,jk) = 0.e0 ! south |
---|
| 135 | ENDIF |
---|
[3] | 136 | |
---|
| 137 | ! ! -------- |
---|
| 138 | ! relative vorticity ! rot |
---|
| 139 | ! ! -------- |
---|
| 140 | ! contravariant velocity (extended for lateral b.c.) |
---|
| 141 | ! inside the model domain |
---|
| 142 | DO jj = 1, jpj |
---|
| 143 | DO ji = 1, jpi |
---|
| 144 | zwu(ji,jj) = e1u(ji,jj) * un(ji,jj,jk) |
---|
| 145 | zwv(ji,jj) = e2v(ji,jj) * vn(ji,jj,jk) |
---|
| 146 | END DO |
---|
| 147 | END DO |
---|
| 148 | |
---|
| 149 | ! East-West boundary conditions |
---|
| 150 | IF( nperio == 1 .OR. nperio == 4 .OR. nperio == 6) THEN |
---|
| 151 | zwv( 0 ,:) = zwv(jpi-2,:) |
---|
| 152 | zwv( -1 ,:) = zwv(jpi-3,:) |
---|
| 153 | zwv(jpi+1,:) = zwv( 3 ,:) |
---|
| 154 | zwv(jpi+2,:) = zwv( 4 ,:) |
---|
| 155 | ELSE |
---|
| 156 | zwv( 0 ,:) = 0.e0 |
---|
| 157 | zwv( -1 ,:) = 0.e0 |
---|
| 158 | zwv(jpi+1,:) = 0.e0 |
---|
| 159 | zwv(jpi+2,:) = 0.e0 |
---|
| 160 | ENDIF |
---|
| 161 | |
---|
| 162 | ! North-South boundary conditions |
---|
| 163 | IF( nperio == 3 .OR. nperio == 4 ) THEN |
---|
| 164 | ! north fold ( Grid defined with a T-point pivot) ORCA 2 degre |
---|
| 165 | zwu(jpi,jpj+1) = 0.e0 |
---|
| 166 | zwu(jpi,jpj+2) = 0.e0 |
---|
| 167 | DO ji = 1, jpi-1 |
---|
| 168 | iju = jpi - ji + 1 |
---|
| 169 | zwu(ji,jpj+1) = - zwu(iju,jpj-3) |
---|
| 170 | zwu(ji,jpj+2) = - zwu(iju,jpj-4) |
---|
| 171 | END DO |
---|
| 172 | ELSEIF( nperio == 5 .OR. nperio == 6 ) THEN |
---|
| 173 | ! north fold ( Grid defined with a F-point pivot) ORCA 0.5 degre\ |
---|
| 174 | zwu(jpi,jpj+1) = 0.e0 |
---|
| 175 | zwu(jpi,jpj+2) = 0.e0 |
---|
| 176 | DO ji = 1, jpi-1 |
---|
| 177 | iju = jpi - ji |
---|
| 178 | zwu(ji,jpj ) = - zwu(iju,jpj-1) |
---|
| 179 | zwu(ji,jpj+1) = - zwu(iju,jpj-2) |
---|
| 180 | zwu(ji,jpj+2) = - zwu(iju,jpj-3) |
---|
| 181 | END DO |
---|
| 182 | DO ji = -1, jpi+2 |
---|
| 183 | ijt = jpi - ji + 1 |
---|
| 184 | zwv(ji,jpj) = - zwv(ijt,jpj-2) |
---|
| 185 | END DO |
---|
| 186 | DO ji = jpi/2+1, jpi+2 |
---|
| 187 | ijt = jpi - ji + 1 |
---|
| 188 | zwv(ji,jpjm1) = - zwv(ijt,jpjm1) |
---|
| 189 | END DO |
---|
| 190 | ELSE |
---|
| 191 | ! closed |
---|
| 192 | zwu(:,jpj+1) = 0.e0 |
---|
| 193 | zwu(:,jpj+2) = 0.e0 |
---|
| 194 | ENDIF |
---|
| 195 | |
---|
| 196 | ! relative vorticity (vertical component of the velocity curl) |
---|
| 197 | DO jj = 1, jpjm1 |
---|
| 198 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
| 199 | rotn(ji,jj,jk) = ( zwv(ji+1,jj ) - zwv(ji,jj) & |
---|
[2528] | 200 | & - zwu(ji ,jj+1) + zwu(ji,jj) ) * fmask(ji,jj,jk) / ( e1f(ji,jj)*e2f(ji,jj) ) |
---|
[3] | 201 | END DO |
---|
| 202 | END DO |
---|
| 203 | |
---|
| 204 | ! second order accurate scheme along straight coast |
---|
| 205 | DO jl = 1, npcoa(1,jk) |
---|
| 206 | ii = nicoa(jl,1,jk) |
---|
| 207 | ij = njcoa(jl,1,jk) |
---|
| 208 | rotn(ii,ij,jk) = 1. / ( e1f(ii,ij) * e2f(ii,ij) ) & |
---|
| 209 | * ( + 4. * zwv(ii+1,ij) - zwv(ii+2,ij) + 0.2 * zwv(ii+3,ij) ) |
---|
| 210 | END DO |
---|
| 211 | DO jl = 1, npcoa(2,jk) |
---|
| 212 | ii = nicoa(jl,2,jk) |
---|
| 213 | ij = njcoa(jl,2,jk) |
---|
| 214 | rotn(ii,ij,jk) = 1./(e1f(ii,ij)*e2f(ii,ij)) & |
---|
| 215 | *(-4.*zwv(ii,ij)+zwv(ii-1,ij)-0.2*zwv(ii-2,ij)) |
---|
| 216 | END DO |
---|
| 217 | DO jl = 1, npcoa(3,jk) |
---|
| 218 | ii = nicoa(jl,3,jk) |
---|
| 219 | ij = njcoa(jl,3,jk) |
---|
| 220 | rotn(ii,ij,jk) = -1. / ( e1f(ii,ij)*e2f(ii,ij) ) & |
---|
| 221 | * ( +4. * zwu(ii,ij+1) - zwu(ii,ij+2) + 0.2 * zwu(ii,ij+3) ) |
---|
| 222 | END DO |
---|
| 223 | DO jl = 1, npcoa(4,jk) |
---|
| 224 | ii = nicoa(jl,4,jk) |
---|
| 225 | ij = njcoa(jl,4,jk) |
---|
| 226 | rotn(ii,ij,jk) = -1. / ( e1f(ii,ij)*e2f(ii,ij) ) & |
---|
| 227 | * ( -4. * zwu(ii,ij) + zwu(ii,ij-1) - 0.2 * zwu(ii,ij-2) ) |
---|
| 228 | END DO |
---|
| 229 | ! ! =============== |
---|
| 230 | END DO ! End of slab |
---|
| 231 | ! ! =============== |
---|
[2528] | 232 | |
---|
[4990] | 233 | IF( ln_rnf ) CALL sbc_rnf_div( hdivn ) ! runoffs (update hdivn field) |
---|
| 234 | IF( ln_divisf .AND. (nn_isf /= 0) ) CALL sbc_isf_div( hdivn ) ! ice shelf (update hdivn field) |
---|
| 235 | IF( nn_cla == 1 ) CALL cla_div ( kt ) ! Cross Land Advection (Update Hor. divergence) |
---|
[3] | 236 | |
---|
| 237 | ! 4. Lateral boundary conditions on hdivn and rotn |
---|
| 238 | ! ---------------------------------=======---====== |
---|
[2528] | 239 | CALL lbc_lnk( hdivn, 'T', 1. ) ; CALL lbc_lnk( rotn , 'F', 1. ) ! lateral boundary cond. (no sign change) |
---|
| 240 | ! |
---|
[6487] | 241 | CALL wrk_dealloc( jpi , jpj+2, zwu ) |
---|
| 242 | CALL wrk_dealloc( jpi+2, jpj , zwv ) |
---|
[3294] | 243 | ! |
---|
| 244 | IF( nn_timing == 1 ) CALL timing_stop('div_cur') |
---|
| 245 | ! |
---|
[3] | 246 | END SUBROUTINE div_cur |
---|
| 247 | |
---|
| 248 | #else |
---|
| 249 | !!---------------------------------------------------------------------- |
---|
| 250 | !! Default option 2nd order centered schemes |
---|
| 251 | !!---------------------------------------------------------------------- |
---|
| 252 | |
---|
| 253 | SUBROUTINE div_cur( kt ) |
---|
| 254 | !!---------------------------------------------------------------------- |
---|
| 255 | !! *** ROUTINE div_cur *** |
---|
| 256 | !! |
---|
| 257 | !! ** Purpose : compute the horizontal divergence and the relative |
---|
| 258 | !! vorticity at before and now time-step |
---|
| 259 | !! |
---|
| 260 | !! ** Method : - Divergence: |
---|
| 261 | !! - save the divergence computed at the previous time-step |
---|
| 262 | !! (note that the Asselin filter has not been applied on hdivb) |
---|
| 263 | !! - compute the now divergence given by : |
---|
| 264 | !! hdivn = 1/(e1t*e2t*e3t) ( di[e2u*e3u un] + dj[e1v*e3v vn] ) |
---|
[2528] | 265 | !! correct hdiv with runoff inflow (div_rnf) and cross land flow (div_cla) |
---|
[3] | 266 | !! - Relavtive Vorticity : |
---|
| 267 | !! - save the curl computed at the previous time-step (rotb = rotn) |
---|
| 268 | !! (note that the Asselin time filter has not been applied to rotb) |
---|
| 269 | !! - compute the now curl in tensorial formalism: |
---|
| 270 | !! rotn = 1/(e1f*e2f) ( di[e2v vn] - dj[e1u un] ) |
---|
| 271 | !! Note: Coastal boundary condition: lateral friction set through |
---|
| 272 | !! the value of fmask along the coast (see dommsk.F90) and shlat |
---|
| 273 | !! (namelist parameter) |
---|
| 274 | !! |
---|
| 275 | !! ** Action : - update hdivb, hdivn, the before & now hor. divergence |
---|
| 276 | !! - update rotb , rotn , the before & now rel. vorticity |
---|
| 277 | !!---------------------------------------------------------------------- |
---|
[2715] | 278 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
[2528] | 279 | ! |
---|
[2715] | 280 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 281 | REAL(wp) :: zraur, zdep ! local scalars |
---|
[3] | 282 | !!---------------------------------------------------------------------- |
---|
[3294] | 283 | ! |
---|
| 284 | IF( nn_timing == 1 ) CALL timing_start('div_cur') |
---|
| 285 | ! |
---|
[3] | 286 | IF( kt == nit000 ) THEN |
---|
| 287 | IF(lwp) WRITE(numout,*) |
---|
| 288 | IF(lwp) WRITE(numout,*) 'div_cur : horizontal velocity divergence and' |
---|
| 289 | IF(lwp) WRITE(numout,*) '~~~~~~~ relative vorticity' |
---|
| 290 | ENDIF |
---|
| 291 | |
---|
| 292 | ! ! =============== |
---|
| 293 | DO jk = 1, jpkm1 ! Horizontal slab |
---|
| 294 | ! ! =============== |
---|
[2528] | 295 | ! |
---|
[3] | 296 | hdivb(:,:,jk) = hdivn(:,:,jk) ! time swap of div arrays |
---|
| 297 | rotb (:,:,jk) = rotn (:,:,jk) ! time swap of rot arrays |
---|
[2528] | 298 | ! |
---|
[3] | 299 | ! ! -------- |
---|
| 300 | ! Horizontal divergence ! div |
---|
| 301 | ! ! -------- |
---|
| 302 | DO jj = 2, jpjm1 |
---|
| 303 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[455] | 304 | hdivn(ji,jj,jk) = & |
---|
[2528] | 305 | ( e2u(ji,jj)*fse3u(ji,jj,jk) * un(ji,jj,jk) - e2u(ji-1,jj)*fse3u(ji-1,jj,jk) * un(ji-1,jj,jk) & |
---|
| 306 | + e1v(ji,jj)*fse3v(ji,jj,jk) * vn(ji,jj,jk) - e1v(ji,jj-1)*fse3v(ji,jj-1,jk) * vn(ji,jj-1,jk) ) & |
---|
[455] | 307 | / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) |
---|
[3] | 308 | END DO |
---|
| 309 | END DO |
---|
| 310 | |
---|
[1792] | 311 | IF( .NOT. AGRIF_Root() ) THEN |
---|
[780] | 312 | IF ((nbondi == 1).OR.(nbondi == 2)) hdivn(nlci-1 , : ,jk) = 0.e0 ! east |
---|
| 313 | IF ((nbondi == -1).OR.(nbondi == 2)) hdivn(2 , : ,jk) = 0.e0 ! west |
---|
| 314 | IF ((nbondj == 1).OR.(nbondj == 2)) hdivn(: ,nlcj-1 ,jk) = 0.e0 ! north |
---|
| 315 | IF ((nbondj == -1).OR.(nbondj == 2)) hdivn(: ,2 ,jk) = 0.e0 ! south |
---|
[1792] | 316 | ENDIF |
---|
[780] | 317 | |
---|
[3] | 318 | ! ! -------- |
---|
| 319 | ! relative vorticity ! rot |
---|
| 320 | ! ! -------- |
---|
| 321 | DO jj = 1, jpjm1 |
---|
| 322 | DO ji = 1, fs_jpim1 ! vector opt. |
---|
| 323 | rotn(ji,jj,jk) = ( e2v(ji+1,jj ) * vn(ji+1,jj ,jk) - e2v(ji,jj) * vn(ji,jj,jk) & |
---|
| 324 | & - e1u(ji ,jj+1) * un(ji ,jj+1,jk) + e1u(ji,jj) * un(ji,jj,jk) ) & |
---|
| 325 | & * fmask(ji,jj,jk) / ( e1f(ji,jj) * e2f(ji,jj) ) |
---|
| 326 | END DO |
---|
| 327 | END DO |
---|
| 328 | ! ! =============== |
---|
| 329 | END DO ! End of slab |
---|
| 330 | ! ! =============== |
---|
[2528] | 331 | |
---|
[4990] | 332 | IF( ln_rnf ) CALL sbc_rnf_div( hdivn ) ! runoffs (update hdivn field) |
---|
| 333 | IF( ln_divisf .AND. (nn_isf .GT. 0) ) CALL sbc_isf_div( hdivn ) ! ice shelf (update hdivn field) |
---|
[2528] | 334 | IF( nn_cla == 1 ) CALL cla_div ( kt ) ! Cross Land Advection (update hdivn field) |
---|
[2715] | 335 | ! |
---|
[2528] | 336 | CALL lbc_lnk( hdivn, 'T', 1. ) ; CALL lbc_lnk( rotn , 'F', 1. ) ! lateral boundary cond. (no sign change) |
---|
| 337 | ! |
---|
[3294] | 338 | IF( nn_timing == 1 ) CALL timing_stop('div_cur') |
---|
| 339 | ! |
---|
[3] | 340 | END SUBROUTINE div_cur |
---|
| 341 | |
---|
| 342 | #endif |
---|
[7649] | 343 | |
---|
| 344 | |
---|
[7845] | 345 | SUBROUTINE dyn_vrt_dia_3d( utend, vtend, id_dia_vor_int, id_dia_vor_mn) |
---|
[7649] | 346 | |
---|
| 347 | !!---------------------------------------------------------------------- |
---|
| 348 | !! |
---|
[7845] | 349 | !! ** Purpose : compute the vertical integrals of utend and vtend, and |
---|
| 350 | !! then pass to dyn_vrt_dia_2d to calculate vorticity |
---|
| 351 | !! tendencies. |
---|
[7649] | 352 | !! |
---|
| 353 | !! ** Action : a) Calculate the vertical integrals of utend & of vtend |
---|
| 354 | !! (u_int & v_int) |
---|
[7845] | 355 | !! b) Call dyn_vrt_dia_2d with vertical integrals |
---|
[7649] | 356 | !! |
---|
| 357 | !!---------------------------------------------------------------------- |
---|
| 358 | REAL :: utend(jpi,jpj,jpk) ! contribution to du/dt |
---|
| 359 | REAL :: vtend(jpi,jpj,jpk) ! contribution to dv/dt |
---|
| 360 | INTEGER :: id_dia_vor_int ! identifier for the vertical integral vorticity diagnostic |
---|
| 361 | INTEGER :: id_dia_vor_mn ! identifier for the vertical mean vorticity diagnostic |
---|
| 362 | ! |
---|
| 363 | !!---------------------------------------------------------------------- |
---|
| 364 | ! |
---|
| 365 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 366 | ! |
---|
| 367 | REAL(wp), POINTER, DIMENSION(:,:) :: u_int ! u vertical integral |
---|
| 368 | REAL(wp), POINTER, DIMENSION(:,:) :: v_int ! v vertical integral |
---|
| 369 | |
---|
| 370 | CALL wrk_alloc(jpi, jpj, u_int) |
---|
| 371 | CALL wrk_alloc(jpi, jpj, v_int) |
---|
| 372 | |
---|
| 373 | u_int(:,:) = 0.0_wp |
---|
| 374 | v_int(:,:) = 0.0_wp |
---|
| 375 | |
---|
| 376 | ! |
---|
| 377 | ! Calculate the vertical integrals of utend & of vtend |
---|
| 378 | ! |
---|
| 379 | |
---|
| 380 | DO jk = 1,jpk |
---|
| 381 | DO jj = 1,jpj |
---|
| 382 | DO ji = 1,jpi |
---|
| 383 | u_int(ji,jj) = u_int(ji,jj) + utend(ji,jj,jk)*fse3u(ji,jj,jk) |
---|
| 384 | v_int(ji,jj) = v_int(ji,jj) + vtend(ji,jj,jk)*fse3v(ji,jj,jk) |
---|
| 385 | END DO |
---|
| 386 | END DO |
---|
| 387 | END DO |
---|
| 388 | |
---|
[7845] | 389 | CALL dyn_vrt_dia_2d(u_int, v_int, id_dia_vor_int, id_dia_vor_mn) |
---|
| 390 | |
---|
| 391 | CALL wrk_dealloc(jpi, jpj, u_int) |
---|
| 392 | CALL wrk_dealloc(jpi, jpj, v_int) |
---|
| 393 | |
---|
| 394 | END SUBROUTINE dyn_vrt_dia_3d |
---|
| 395 | |
---|
| 396 | |
---|
| 397 | SUBROUTINE dyn_vrt_dia_2d( u_int, v_int, id_dia_vor_int, id_dia_vor_mn) |
---|
| 398 | |
---|
| 399 | !!---------------------------------------------------------------------- |
---|
| 400 | !! |
---|
| 401 | !! ** Purpose : compute the integral and mean vorticity tendencies. |
---|
| 402 | !! |
---|
| 403 | !! ** Action : a) Calculate the vorticity tendencies for the vertical |
---|
| 404 | !! integrals. |
---|
| 405 | !! b) Calculate the vertical means, u_mn, v_mn from u_int |
---|
| 406 | !! and v_int by dividing by the depth |
---|
| 407 | !! c) Calculate the vorticity tendencies for the vertical |
---|
| 408 | !! means |
---|
| 409 | !! d) Call iom_put for the vertical integral vorticity |
---|
| 410 | !! tendencies (using id_dia_vor_int) |
---|
| 411 | !! e) Call iom_put for the vertical mean vorticity |
---|
| 412 | !! tendencies (using id_dia_vor_mn) |
---|
| 413 | !! |
---|
| 414 | !!---------------------------------------------------------------------- |
---|
| 415 | REAL :: u_int(jpi,jpj) ! u vertical integral |
---|
| 416 | REAL :: v_int(jpi,jpj) ! v vertical integral |
---|
| 417 | INTEGER :: id_dia_vor_int ! identifier for the vertical integral vorticity diagnostic |
---|
| 418 | INTEGER :: id_dia_vor_mn ! identifier for the vertical mean vorticity diagnostic |
---|
[7649] | 419 | ! |
---|
[7845] | 420 | !!---------------------------------------------------------------------- |
---|
| 421 | ! |
---|
| 422 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 423 | INTEGER :: ikbu, ikbv ! dummy loop indices |
---|
| 424 | ! |
---|
| 425 | REAL(wp), POINTER, DIMENSION(:,:) :: u_mn ! u vertical means |
---|
| 426 | REAL(wp), POINTER, DIMENSION(:,:) :: v_mn ! u vertical means |
---|
| 427 | REAL(wp), POINTER, DIMENSION(:,:) :: vor_int ! vort trend of vert integrals |
---|
| 428 | REAL(wp), POINTER, DIMENSION(:,:) :: vor_mn ! vort trend of vert means |
---|
| 429 | |
---|
| 430 | CALL wrk_alloc(jpi, jpj, u_mn) |
---|
| 431 | CALL wrk_alloc(jpi, jpj, v_mn) |
---|
| 432 | CALL wrk_alloc(jpi, jpj, vor_int) |
---|
| 433 | CALL wrk_alloc(jpi, jpj, vor_mn) |
---|
| 434 | |
---|
[7900] | 435 | CALL lbc_lnk( u_int, 'U', 1. ) |
---|
| 436 | CALL lbc_lnk( v_int, 'V', 1. ) |
---|
| 437 | |
---|
[7845] | 438 | ! |
---|
[7649] | 439 | ! Calculate the vorticity tendencies for the vertical integrals. |
---|
| 440 | ! 1/e1e2 * ((e2*d(vtend)/dx) - (e1*d(utend)/dy)) |
---|
| 441 | ! |
---|
| 442 | |
---|
| 443 | DO jj = 1,jpjm1 |
---|
| 444 | DO ji = 1,jpim1 |
---|
| 445 | vor_int(ji,jj) = ( v_int(ji+1,jj) * e2v(ji+1,jj) & |
---|
| 446 | & - v_int(ji,jj) * e2v(ji,jj) & |
---|
| 447 | & + u_int(ji,jj) * e1u(ji,jj) & |
---|
| 448 | & - u_int(ji,jj+1) * e1u(ji,jj+1) ) & |
---|
| 449 | & / ( e1f(ji,jj) * e2f(ji,jj) ) |
---|
| 450 | END DO |
---|
| 451 | END DO |
---|
| 452 | |
---|
| 453 | |
---|
| 454 | ! |
---|
| 455 | ! Calculate the vertical means, u_mn, v_mn from u_int & v_int by dividing |
---|
| 456 | ! by the depth |
---|
| 457 | ! mbku & mbkv - vertical index of the bottom last U- & W- ocean level |
---|
| 458 | ! |
---|
| 459 | |
---|
| 460 | DO jj = 1, jpj |
---|
| 461 | DO ji = 1, jpi |
---|
| 462 | ikbu = mbku(ji,jj) |
---|
| 463 | ikbv = mbkv(ji,jj) |
---|
| 464 | |
---|
| 465 | IF (ikbu .ne. 0.0_wp) THEN ! Don't divide by 0! |
---|
| 466 | u_mn(ji,jj) = u_int(ji,jj) / ikbu |
---|
| 467 | ELSE |
---|
| 468 | u_mn(ji,jj) = 0.0_wp |
---|
| 469 | END IF |
---|
| 470 | |
---|
| 471 | IF (ikbv .ne. 0.0_wp) THEN ! Don't divide by 0! |
---|
| 472 | v_mn(ji,jj) = v_int(ji,jj) / ikbv |
---|
| 473 | ELSE |
---|
| 474 | v_mn(ji,jj) = 0.0_wp |
---|
| 475 | END IF |
---|
| 476 | END DO |
---|
| 477 | END DO |
---|
| 478 | |
---|
| 479 | ! |
---|
| 480 | ! Calculate the vorticity tendencies for the vertical means |
---|
| 481 | ! 1/e1e2 * ((e2*d(v_mn)/dx) - (e1*d(u_mn)/dy)) |
---|
| 482 | ! |
---|
| 483 | |
---|
| 484 | DO jj = 1,jpjm1 |
---|
| 485 | DO ji = 1,jpim1 |
---|
| 486 | vor_mn(ji,jj) = ( v_mn(ji+1,jj) * e2v(ji+1,jj) & |
---|
| 487 | & - v_mn(ji,jj) * e2v(ji,jj) & |
---|
| 488 | & + u_mn(ji,jj) * e1u(ji,jj) & |
---|
| 489 | & - u_mn(ji,jj+1) * e1u(ji,jj+1) ) & |
---|
| 490 | & / ( e1f(ji,jj) * e2f(ji,jj) ) |
---|
| 491 | END DO |
---|
| 492 | END DO |
---|
| 493 | |
---|
| 494 | |
---|
| 495 | ! Call iom_put for the vertical integral vorticity tendencies |
---|
| 496 | IF (id_dia_vor_int == 1) THEN |
---|
| 497 | CALL iom_put( "dia_vor_int", vor_int(:,:)) |
---|
| 498 | ENDIF |
---|
| 499 | |
---|
| 500 | ! Call iom_put for the vertical mean vorticity tendencies |
---|
| 501 | IF (id_dia_vor_int == 1) THEN |
---|
| 502 | CALL iom_put( "dia_vor_mn", vor_mn(:,:)) |
---|
| 503 | ENDIF |
---|
| 504 | |
---|
| 505 | CALL wrk_dealloc(jpi, jpj, u_mn) |
---|
| 506 | CALL wrk_dealloc(jpi, jpj, v_mn) |
---|
| 507 | CALL wrk_dealloc(jpi, jpj, vor_int) |
---|
| 508 | CALL wrk_dealloc(jpi, jpj, vor_mn) |
---|
| 509 | |
---|
[7845] | 510 | END SUBROUTINE dyn_vrt_dia_2d |
---|
[7649] | 511 | |
---|
| 512 | |
---|
[3] | 513 | !!====================================================================== |
---|
| 514 | END MODULE divcur |
---|