| 1 | #include "cppdefs.h"
|
|---|
| 2 | MODULE set_avg_mod
|
|---|
| 3 | #if defined AVERAGES && (!defined ADJOINT && defined NONLINEAR)
|
|---|
| 4 | !
|
|---|
| 5 | !svn $Id: set_avg.F 42 2007-05-04 20:58:11Z arango $
|
|---|
| 6 | !================================================== Hernan G. Arango ===
|
|---|
| 7 | ! Copyright (c) 2002-2007 The ROMS/TOMS Group !
|
|---|
| 8 | ! Licensed under a MIT/X style license !
|
|---|
| 9 | ! See License_ROMS.txt !
|
|---|
| 10 | !=======================================================================
|
|---|
| 11 | ! !
|
|---|
| 12 | ! This subroutine accumulates and computes output time-averaged !
|
|---|
| 13 | ! fields. Due to synchronization, the time-averaged fields are !
|
|---|
| 14 | ! computed in delayed mode. All averages are accumulated at the !
|
|---|
| 15 | ! beggining of the next time-step. !
|
|---|
| 16 | ! !
|
|---|
| 17 | !=======================================================================
|
|---|
| 18 | !
|
|---|
| 19 | implicit none
|
|---|
| 20 |
|
|---|
| 21 | PRIVATE
|
|---|
| 22 | PUBLIC :: set_avg
|
|---|
| 23 |
|
|---|
| 24 | CONTAINS
|
|---|
| 25 | !
|
|---|
| 26 | !***********************************************************************
|
|---|
| 27 | SUBROUTINE set_avg (ng, tile)
|
|---|
| 28 | !***********************************************************************
|
|---|
| 29 | !
|
|---|
| 30 | USE mod_param
|
|---|
| 31 | USE mod_average
|
|---|
| 32 | USE mod_forces
|
|---|
| 33 | # ifdef SOLVE3D
|
|---|
| 34 | USE mod_grid
|
|---|
| 35 | USE mod_mixing
|
|---|
| 36 | # endif
|
|---|
| 37 | USE mod_ocean
|
|---|
| 38 | USE mod_stepping
|
|---|
| 39 | !
|
|---|
| 40 | ! Imported variable declarations.
|
|---|
| 41 | !
|
|---|
| 42 | integer, intent(in) :: ng, tile
|
|---|
| 43 | !
|
|---|
| 44 | ! Local variable declarations.
|
|---|
| 45 | !
|
|---|
| 46 | # include "tile.h"
|
|---|
| 47 | !
|
|---|
| 48 | # ifdef PROFILE
|
|---|
| 49 | CALL wclock_on (ng, iNLM, 5)
|
|---|
| 50 | # endif
|
|---|
| 51 | CALL set_avg_tile (ng, Istr, Iend, Jstr, Jend, &
|
|---|
| 52 | & LBi, UBi, LBj, UBj, &
|
|---|
| 53 | & KOUT, &
|
|---|
| 54 | # ifdef SOLVE3D
|
|---|
| 55 | & NOUT, &
|
|---|
| 56 | & GRID(ng) % pm, &
|
|---|
| 57 | & GRID(ng) % pn, &
|
|---|
| 58 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 59 | & GRID(ng) % Huon, &
|
|---|
| 60 | & GRID(ng) % Hvom, &
|
|---|
| 61 | # endif
|
|---|
| 62 | & OCEAN(ng) % u, &
|
|---|
| 63 | & OCEAN(ng) % v, &
|
|---|
| 64 | & OCEAN(ng) % W, &
|
|---|
| 65 | & OCEAN(ng) % wvel, &
|
|---|
| 66 | & OCEAN(ng) % t, &
|
|---|
| 67 | & OCEAN(ng) % rho, &
|
|---|
| 68 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 69 | & OCEAN(ng) % u_stokes, &
|
|---|
| 70 | & OCEAN(ng) % v_stokes, &
|
|---|
| 71 | & MIXING(ng) % rustr3d, &
|
|---|
| 72 | & MIXING(ng) % rvstr3d, &
|
|---|
| 73 | & MIXING(ng) % Sxx, &
|
|---|
| 74 | & MIXING(ng) % Sxy, &
|
|---|
| 75 | & MIXING(ng) % Syy, &
|
|---|
| 76 | & MIXING(ng) % Szx, &
|
|---|
| 77 | & MIXING(ng) % Szy, &
|
|---|
| 78 | # endif
|
|---|
| 79 | # ifdef LMD_SKPP
|
|---|
| 80 | & MIXING(ng) % hsbl, &
|
|---|
| 81 | # endif
|
|---|
| 82 | # ifdef LMD_BKPP
|
|---|
| 83 | & MIXING(ng) % hbbl, &
|
|---|
| 84 | # endif
|
|---|
| 85 | # ifdef AVERAGES_AKV
|
|---|
| 86 | & MIXING(ng) % Akv, &
|
|---|
| 87 | # endif
|
|---|
| 88 | # if defined AVERAGES_AKT || defined AVERAGES_AKS
|
|---|
| 89 | & MIXING(ng) % Akt, &
|
|---|
| 90 | # endif
|
|---|
| 91 | # ifdef AVERAGES_FLUXES
|
|---|
| 92 | & FORCES(ng) % stflx, &
|
|---|
| 93 | # ifdef BULK_FLUXES
|
|---|
| 94 | & FORCES(ng) % lhflx, &
|
|---|
| 95 | & FORCES(ng) % shflx, &
|
|---|
| 96 | & FORCES(ng) % lrflx, &
|
|---|
| 97 | # ifdef EMINUSP
|
|---|
| 98 | & FORCES(ng) % evap, &
|
|---|
| 99 | & FORCES(ng) % rain, &
|
|---|
| 100 | # endif
|
|---|
| 101 | # endif
|
|---|
| 102 | # ifdef SHORTWAVE
|
|---|
| 103 | & FORCES(ng) % srflx, &
|
|---|
| 104 | # endif
|
|---|
| 105 | # endif
|
|---|
| 106 | # endif
|
|---|
| 107 | # ifdef AVERAGES_FLUXES
|
|---|
| 108 | & FORCES(ng) % sustr, &
|
|---|
| 109 | & FORCES(ng) % svstr, &
|
|---|
| 110 | # endif
|
|---|
| 111 | & OCEAN(ng) % ubar, &
|
|---|
| 112 | & OCEAN(ng) % vbar, &
|
|---|
| 113 | & OCEAN(ng) % zeta, &
|
|---|
| 114 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 115 | & OCEAN(ng) % ubar_stokes, &
|
|---|
| 116 | & OCEAN(ng) % vbar_stokes, &
|
|---|
| 117 | & MIXING(ng) % rustr2d, &
|
|---|
| 118 | & MIXING(ng) % rvstr2d, &
|
|---|
| 119 | & MIXING(ng) % Sxx_bar, &
|
|---|
| 120 | & MIXING(ng) % Sxy_bar, &
|
|---|
| 121 | & MIXING(ng) % Syy_bar, &
|
|---|
| 122 | # endif
|
|---|
| 123 | # ifdef SOLVE3D
|
|---|
| 124 | & AVERAGE(ng) % avgu3d, &
|
|---|
| 125 | & AVERAGE(ng) % avgv3d, &
|
|---|
| 126 | & AVERAGE(ng) % avgw3d, &
|
|---|
| 127 | & AVERAGE(ng) % avgwvel, &
|
|---|
| 128 | & AVERAGE(ng) % avgt, &
|
|---|
| 129 | & AVERAGE(ng) % avgrho, &
|
|---|
| 130 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 131 | & AVERAGE(ng) % avgu3Sd, &
|
|---|
| 132 | & AVERAGE(ng) % avgv3Sd, &
|
|---|
| 133 | & AVERAGE(ng) % avgu3RS, &
|
|---|
| 134 | & AVERAGE(ng) % avgv3RS, &
|
|---|
| 135 | & AVERAGE(ng) % avgSxx3d, &
|
|---|
| 136 | & AVERAGE(ng) % avgSxy3d, &
|
|---|
| 137 | & AVERAGE(ng) % avgSyy3d, &
|
|---|
| 138 | & AVERAGE(ng) % avgSzx3d, &
|
|---|
| 139 | & AVERAGE(ng) % avgSzy3d, &
|
|---|
| 140 | # endif
|
|---|
| 141 | # ifdef LMD_SKPP
|
|---|
| 142 | & AVERAGE(ng) % avghsbl, &
|
|---|
| 143 | # endif
|
|---|
| 144 | # ifdef LMD_BKPP
|
|---|
| 145 | & AVERAGE(ng) % avghbbl, &
|
|---|
| 146 | # endif
|
|---|
| 147 | # ifdef AVERAGES_AKV
|
|---|
| 148 | & AVERAGE(ng) % avgAKv, &
|
|---|
| 149 | # endif
|
|---|
| 150 | # ifdef AVERAGES_AKT
|
|---|
| 151 | & AVERAGE(ng) % avgAKt, &
|
|---|
| 152 | # endif
|
|---|
| 153 | # ifdef AVERAGES_AKS
|
|---|
| 154 | & AVERAGE(ng) % avgAKs, &
|
|---|
| 155 | # endif
|
|---|
| 156 | # ifdef AVERAGES_FLUXES
|
|---|
| 157 | & AVERAGE(ng) % avgstf, &
|
|---|
| 158 | & AVERAGE(ng) % avgswf, &
|
|---|
| 159 | # ifdef BULK_FLUXES
|
|---|
| 160 | & AVERAGE(ng) % avglhf, &
|
|---|
| 161 | & AVERAGE(ng) % avgshf, &
|
|---|
| 162 | & AVERAGE(ng) % avglrf, &
|
|---|
| 163 | # ifdef EMINUSP
|
|---|
| 164 | & AVERAGE(ng) % avgevap, &
|
|---|
| 165 | & AVERAGE(ng) % avgrain, &
|
|---|
| 166 | # endif
|
|---|
| 167 | # endif
|
|---|
| 168 | # ifdef SHORTWAVE
|
|---|
| 169 | & AVERAGE(ng) % avgsrf, &
|
|---|
| 170 | # endif
|
|---|
| 171 | # endif
|
|---|
| 172 | # endif
|
|---|
| 173 | # ifdef AVERAGES_FLUXES
|
|---|
| 174 | & AVERAGE(ng) % avgsus, &
|
|---|
| 175 | & AVERAGE(ng) % avgsvs, &
|
|---|
| 176 | # endif
|
|---|
| 177 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 178 | # ifdef SOLVE3D
|
|---|
| 179 | & AVERAGE(ng) % avgHuon, &
|
|---|
| 180 | & AVERAGE(ng) % avgHvom, &
|
|---|
| 181 | & AVERAGE(ng) % avgUU, &
|
|---|
| 182 | & AVERAGE(ng) % avgUV, &
|
|---|
| 183 | & AVERAGE(ng) % avgVV, &
|
|---|
| 184 | & AVERAGE(ng) % avgHuonT, &
|
|---|
| 185 | & AVERAGE(ng) % avgHvomT, &
|
|---|
| 186 | & AVERAGE(ng) % avgUT, &
|
|---|
| 187 | & AVERAGE(ng) % avgVT, &
|
|---|
| 188 | & AVERAGE(ng) % avgTT, &
|
|---|
| 189 | # endif
|
|---|
| 190 | & AVERAGE(ng) % avgU2, &
|
|---|
| 191 | & AVERAGE(ng) % avgV2, &
|
|---|
| 192 | & AVERAGE(ng) % avgZZ, &
|
|---|
| 193 | # endif
|
|---|
| 194 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 195 | & AVERAGE(ng) % avgu2Sd, &
|
|---|
| 196 | & AVERAGE(ng) % avgv2Sd, &
|
|---|
| 197 | & AVERAGE(ng) % avgu2RS, &
|
|---|
| 198 | & AVERAGE(ng) % avgv2RS, &
|
|---|
| 199 | & AVERAGE(ng) % avgSxx2d, &
|
|---|
| 200 | & AVERAGE(ng) % avgSxy2d, &
|
|---|
| 201 | & AVERAGE(ng) % avgSyy2d, &
|
|---|
| 202 | # endif
|
|---|
| 203 | & AVERAGE(ng) % avgu2d, &
|
|---|
| 204 | & AVERAGE(ng) % avgv2d, &
|
|---|
| 205 | & AVERAGE(ng) % avgzeta)
|
|---|
| 206 |
|
|---|
| 207 | # ifdef PROFILE
|
|---|
| 208 | CALL wclock_off (ng, iNLM, 5)
|
|---|
| 209 | # endif
|
|---|
| 210 | RETURN
|
|---|
| 211 | END SUBROUTINE set_avg
|
|---|
| 212 | !
|
|---|
| 213 | !***********************************************************************
|
|---|
| 214 | SUBROUTINE set_avg_tile (ng, Istr, Iend, Jstr, Jend, &
|
|---|
| 215 | & LBi, UBi, LBj, UBj, &
|
|---|
| 216 | & Kout, &
|
|---|
| 217 | # ifdef SOLVE3D
|
|---|
| 218 | & Nout, &
|
|---|
| 219 | & pm, pn, &
|
|---|
| 220 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 221 | & Huon, Hvom, &
|
|---|
| 222 | # endif
|
|---|
| 223 | & u, v, &
|
|---|
| 224 | & W, wvel, t, rho, &
|
|---|
| 225 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 226 | & u_stokes, v_stokes, &
|
|---|
| 227 | & rustr3d, rvstr3d, &
|
|---|
| 228 | & Sxx, Sxy, Syy, Szx, Szy, &
|
|---|
| 229 | # endif
|
|---|
| 230 | # ifdef LMD_SKPP
|
|---|
| 231 | & hsbl, &
|
|---|
| 232 | # endif
|
|---|
| 233 | # ifdef LMD_BKPP
|
|---|
| 234 | & hbbl, &
|
|---|
| 235 | # endif
|
|---|
| 236 | # ifdef AVERAGES_AKV
|
|---|
| 237 | & Akv, &
|
|---|
| 238 | # endif
|
|---|
| 239 | # if defined AVERAGES_AKT || defined AVERAGES_AKS
|
|---|
| 240 | & Akt, &
|
|---|
| 241 | # endif
|
|---|
| 242 | # ifdef AVERAGES_FLUXES
|
|---|
| 243 | & stflx, &
|
|---|
| 244 | # ifdef BULK_FLUXES
|
|---|
| 245 | & lhflx, shflx, lrflx, &
|
|---|
| 246 | # ifdef EMINUSP
|
|---|
| 247 | & evap, rain, &
|
|---|
| 248 | # endif
|
|---|
| 249 | # endif
|
|---|
| 250 | # ifdef SHORTWAVE
|
|---|
| 251 | & srflx, &
|
|---|
| 252 | # endif
|
|---|
| 253 | # endif
|
|---|
| 254 | # endif
|
|---|
| 255 | # ifdef AVERAGES_FLUXES
|
|---|
| 256 | & sustr, svstr, &
|
|---|
| 257 | # endif
|
|---|
| 258 | & ubar, vbar, &
|
|---|
| 259 | & zeta, &
|
|---|
| 260 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 261 | & ubar_stokes, vbar_stokes, &
|
|---|
| 262 | & rustr2d, rvstr2d, &
|
|---|
| 263 | & Sxx_bar, Sxy_bar, Syy_bar, &
|
|---|
| 264 | # endif
|
|---|
| 265 | # ifdef SOLVE3D
|
|---|
| 266 | & avgu3d, avgv3d, &
|
|---|
| 267 | & avgw3d, avgwvel, &
|
|---|
| 268 | & avgt, avgrho, &
|
|---|
| 269 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 270 | & avgu3Sd, avgv3Sd, &
|
|---|
| 271 | & avgu3RS, avgv3RS, &
|
|---|
| 272 | & avgSxx3d, avgSxy3d, avgSyy3d, &
|
|---|
| 273 | & avgSzx3d, avgSzy3d, &
|
|---|
| 274 | # endif
|
|---|
| 275 | # ifdef LMD_SKPP
|
|---|
| 276 | & avghsbl, &
|
|---|
| 277 | # endif
|
|---|
| 278 | # ifdef LMD_BKPP
|
|---|
| 279 | & avghbbl, &
|
|---|
| 280 | # endif
|
|---|
| 281 | # ifdef AVERAGES_AKV
|
|---|
| 282 | & avgAKv, &
|
|---|
| 283 | # endif
|
|---|
| 284 | # ifdef AVERAGES_AKT
|
|---|
| 285 | & avgAKt, &
|
|---|
| 286 | # endif
|
|---|
| 287 | # ifdef AVERAGES_AKS
|
|---|
| 288 | & avgAKs, &
|
|---|
| 289 | # endif
|
|---|
| 290 | # ifdef AVERAGES_FLUXES
|
|---|
| 291 | & avgstf, avgswf, &
|
|---|
| 292 | # ifdef BULK_FLUXES
|
|---|
| 293 | & avglhf, avgshf, avglrf, &
|
|---|
| 294 | # ifdef EMINUSP
|
|---|
| 295 | & avgevap, avgrain, &
|
|---|
| 296 | # endif
|
|---|
| 297 | # endif
|
|---|
| 298 | # ifdef SHORTWAVE
|
|---|
| 299 | & avgsrf, &
|
|---|
| 300 | # endif
|
|---|
| 301 | # endif
|
|---|
| 302 | # endif
|
|---|
| 303 | # ifdef AVERAGES_FLUXES
|
|---|
| 304 | & avgsus, avgsvs, &
|
|---|
| 305 | # endif
|
|---|
| 306 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 307 | # ifdef SOLVE3D
|
|---|
| 308 | & avgHuon, avgHvom, &
|
|---|
| 309 | & avgUU, avgUV, avgVV, &
|
|---|
| 310 | & avgHuonT, avgHvomT, &
|
|---|
| 311 | & avgUT, avgVT, avgTT, &
|
|---|
| 312 | # endif
|
|---|
| 313 | & avgU2, avgV2, avgZZ, &
|
|---|
| 314 | # endif
|
|---|
| 315 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 316 | & avgu2Sd, avgv2Sd, &
|
|---|
| 317 | & avgu2RS, avgv2RS, &
|
|---|
| 318 | & avgSxx2d, avgSxy2d, avgSyy2d, &
|
|---|
| 319 | # endif
|
|---|
| 320 | & avgu2d, avgv2d, &
|
|---|
| 321 | & avgzeta)
|
|---|
| 322 | !***********************************************************************
|
|---|
| 323 | !
|
|---|
| 324 | USE mod_param
|
|---|
| 325 | USE mod_scalars
|
|---|
| 326 | !
|
|---|
| 327 | implicit none
|
|---|
| 328 | !
|
|---|
| 329 | ! Imported variable declarations.
|
|---|
| 330 | !
|
|---|
| 331 | integer, intent(in) :: ng, Iend, Istr, Jend, Jstr
|
|---|
| 332 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
|---|
| 333 | integer, intent(in) :: Kout
|
|---|
| 334 | # ifdef SOLVE3D
|
|---|
| 335 | integer, intent(in) :: Nout
|
|---|
| 336 | # endif
|
|---|
| 337 | !
|
|---|
| 338 | # ifdef ASSUMED_SHAPE
|
|---|
| 339 | # ifdef SOLVE3D
|
|---|
| 340 | real(r8), intent(in) :: pm(LBi:,LBj:)
|
|---|
| 341 | real(r8), intent(in) :: pn(LBi:,LBj:)
|
|---|
| 342 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 343 | real(r8), intent(in) :: Huon(LBi:,LBj:,:)
|
|---|
| 344 | real(r8), intent(in) :: Hvom(LBi:,LBj:,:)
|
|---|
| 345 | # endif
|
|---|
| 346 | real(r8), intent(in) :: u(LBi:,LBj:,:,:)
|
|---|
| 347 | real(r8), intent(in) :: v(LBi:,LBj:,:,:)
|
|---|
| 348 | real(r8), intent(in) :: W(LBi:,LBj:,0:)
|
|---|
| 349 | real(r8), intent(in) :: wvel(LBi:,LBj:,0:)
|
|---|
| 350 | real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
|
|---|
| 351 | real(r8), intent(in) :: rho(LBi:,LBj:,:)
|
|---|
| 352 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 353 | real(r8), intent(in) :: u_stokes(LBi:,LBj:,:)
|
|---|
| 354 | real(r8), intent(in) :: v_stokes(LBi:,LBj:,:)
|
|---|
| 355 | real(r8), intent(in) :: rustr3d(LBi:,LBj:,:)
|
|---|
| 356 | real(r8), intent(in) :: rvstr3d(LBi:,LBj:,:)
|
|---|
| 357 | real(r8), intent(in) :: Sxx(LBi:,LBj:,:)
|
|---|
| 358 | real(r8), intent(in) :: Sxy(LBi:,LBj:,:)
|
|---|
| 359 | real(r8), intent(in) :: Syy(LBi:,LBj:,:)
|
|---|
| 360 | real(r8), intent(in) :: Szx(LBi:,LBj:,:)
|
|---|
| 361 | real(r8), intent(in) :: Szy(LBi:,LBj:,:)
|
|---|
| 362 | # endif
|
|---|
| 363 | # ifdef LMD_SKPP
|
|---|
| 364 | real(r8), intent(in) :: hsbl(LBi:,LBj:)
|
|---|
| 365 | # endif
|
|---|
| 366 | # ifdef LMD_BKPP
|
|---|
| 367 | real(r8), intent(in) :: hbbl(LBi:,LBj:)
|
|---|
| 368 | # endif
|
|---|
| 369 | # ifdef AVERAGES_AKV
|
|---|
| 370 | real(r8), intent(in) :: Akv(LBi:,LBj:,0:)
|
|---|
| 371 | # endif
|
|---|
| 372 | # if defined AVERAGES_AKT || defined AVERAGES_AKS
|
|---|
| 373 | real(r8), intent(in) :: Akt(LBi:,LBj:,0:,:)
|
|---|
| 374 | # endif
|
|---|
| 375 | # ifdef AVERAGES_FLUXES
|
|---|
| 376 | real(r8), intent(in) :: stflx(LBi:,LBj:,:)
|
|---|
| 377 | # ifdef BULK_FLUXES
|
|---|
| 378 | real(r8), intent(in) :: lhflx(LBi:,LBj:)
|
|---|
| 379 | real(r8), intent(in) :: shflx(LBi:,LBj:)
|
|---|
| 380 | real(r8), intent(in) :: lrflx(LBi:,LBj:)
|
|---|
| 381 | # ifdef EMINUSP
|
|---|
| 382 | real(r8), intent(in) :: evap(LBi:,LBj:)
|
|---|
| 383 | real(r8), intent(in) :: rain(LBi:,LBj:)
|
|---|
| 384 | # endif
|
|---|
| 385 | # endif
|
|---|
| 386 | # ifdef SHORTWAVE
|
|---|
| 387 | real(r8), intent(in) :: srflx(LBi:,LBj:)
|
|---|
| 388 | # endif
|
|---|
| 389 | # endif
|
|---|
| 390 | # endif
|
|---|
| 391 | # ifdef AVERAGES_FLUXES
|
|---|
| 392 | real(r8), intent(in) :: sustr(LBi:,LBj:)
|
|---|
| 393 | real(r8), intent(in) :: svstr(LBi:,LBj:)
|
|---|
| 394 | # endif
|
|---|
| 395 | real(r8), intent(in) :: ubar(LBi:,LBj:,:)
|
|---|
| 396 | real(r8), intent(in) :: vbar(LBi:,LBj:,:)
|
|---|
| 397 | real(r8), intent(in) :: zeta(LBi:,LBj:,:)
|
|---|
| 398 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 399 | real(r8), intent(in) :: ubar_stokes(LBi:,LBj:)
|
|---|
| 400 | real(r8), intent(in) :: vbar_stokes(LBi:,LBj:)
|
|---|
| 401 | real(r8), intent(in) :: rustr2d(LBi:,LBj:)
|
|---|
| 402 | real(r8), intent(in) :: rvstr2d(LBi:,LBj:)
|
|---|
| 403 | real(r8), intent(in) :: Sxx_bar(LBi:,LBj:)
|
|---|
| 404 | real(r8), intent(in) :: Sxy_bar(LBi:,LBj:)
|
|---|
| 405 | real(r8), intent(in) :: Syy_bar(LBi:,LBj:)
|
|---|
| 406 | # endif
|
|---|
| 407 | # ifdef SOLVE3D
|
|---|
| 408 | real(r8), intent(inout) :: avgu3d(LBi:,LBj:,:)
|
|---|
| 409 | real(r8), intent(inout) :: avgv3d(LBi:,LBj:,:)
|
|---|
| 410 | real(r8), intent(inout) :: avgw3d(LBi:,LBj:,0:)
|
|---|
| 411 | real(r8), intent(inout) :: avgwvel(LBi:,LBj:,0:)
|
|---|
| 412 | real(r8), intent(inout) :: avgt(LBi:,LBj:,:,:)
|
|---|
| 413 | real(r8), intent(inout) :: avgrho(LBi:,LBj:,:)
|
|---|
| 414 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 415 | real(r8), intent(inout) :: avgu3Sd(LBi:,LBj:,:)
|
|---|
| 416 | real(r8), intent(inout) :: avgv3Sd(LBi:,LBj:,:)
|
|---|
| 417 | real(r8), intent(inout) :: avgu3RS(LBi:,LBj:,:)
|
|---|
| 418 | real(r8), intent(inout) :: avgv3RS(LBi:,LBj:,:)
|
|---|
| 419 | real(r8), intent(inout) :: avgSxx3d(LBi:,LBj:,:)
|
|---|
| 420 | real(r8), intent(inout) :: avgSxy3d(LBi:,LBj:,:)
|
|---|
| 421 | real(r8), intent(inout) :: avgSyy3d(LBi:,LBj:,:)
|
|---|
| 422 | real(r8), intent(inout) :: avgSzx3d(LBi:,LBj:,:)
|
|---|
| 423 | real(r8), intent(inout) :: avgSzy3d(LBi:,LBj:,:)
|
|---|
| 424 | # endif
|
|---|
| 425 | # ifdef LMD_SKPP
|
|---|
| 426 | real(r8), intent(inout) :: avghsbl(LBi:,LBj:)
|
|---|
| 427 | # endif
|
|---|
| 428 | # ifdef LMD_BKPP
|
|---|
| 429 | real(r8), intent(inout) :: avghbbl(LBi:,LBj:)
|
|---|
| 430 | # endif
|
|---|
| 431 | # ifdef AVERAGES_AKV
|
|---|
| 432 | real(r8), intent(inout) :: avgAKv(LBi:,LBj:,0:)
|
|---|
| 433 | # endif
|
|---|
| 434 | # ifdef AVERAGES_AKT
|
|---|
| 435 | real(r8), intent(inout) :: avgAKt(LBi:,LBj:,0:)
|
|---|
| 436 | # endif
|
|---|
| 437 | # ifdef AVERAGES_AKS
|
|---|
| 438 | real(r8), intent(inout) :: avgAKs(LBi:,LBj:,0:)
|
|---|
| 439 | # endif
|
|---|
| 440 | # ifdef AVERAGES_FLUXES
|
|---|
| 441 | real(r8), intent(inout) :: avgstf(LBi:,LBj:)
|
|---|
| 442 | real(r8), intent(inout) :: avgswf(LBi:,LBj:)
|
|---|
| 443 | # ifdef BULK_FLUXES
|
|---|
| 444 | real(r8), intent(inout) :: avglhf(LBi:,LBj:)
|
|---|
| 445 | real(r8), intent(inout) :: avgshf(LBi:,LBj:)
|
|---|
| 446 | real(r8), intent(inout) :: avglrf(LBi:,LBj:)
|
|---|
| 447 | # ifdef EMINUSP
|
|---|
| 448 | real(r8), intent(inout) :: avgevap(LBi:,LBj:)
|
|---|
| 449 | real(r8), intent(inout) :: avgrain(LBi:,LBj:)
|
|---|
| 450 | # endif
|
|---|
| 451 | # endif
|
|---|
| 452 | # ifdef SHORTWAVE
|
|---|
| 453 | real(r8), intent(inout) :: avgsrf(LBi:,LBj:)
|
|---|
| 454 | # endif
|
|---|
| 455 | # endif
|
|---|
| 456 | # endif
|
|---|
| 457 | # ifdef AVERAGES_FLUXES
|
|---|
| 458 | real(r8), intent(inout) :: avgsus(LBi:,LBj:)
|
|---|
| 459 | real(r8), intent(inout) :: avgsvs(LBi:,LBj:)
|
|---|
| 460 | # endif
|
|---|
| 461 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 462 | # ifdef SOLVE3D
|
|---|
| 463 | real(r8), intent(inout) :: avgHuon(LBi:,LBj:,:)
|
|---|
| 464 | real(r8), intent(inout) :: avgHvom(LBi:,LBj:,:)
|
|---|
| 465 | real(r8), intent(inout) :: avgUU(LBi:,LBj:,:)
|
|---|
| 466 | real(r8), intent(inout) :: avgUV(LBi:,LBj:,:)
|
|---|
| 467 | real(r8), intent(inout) :: avgVV(LBi:,LBj:,:)
|
|---|
| 468 | real(r8), intent(inout) :: avgHuonT(LBi:,LBj:,:,:)
|
|---|
| 469 | real(r8), intent(inout) :: avgHvomT(LBi:,LBj:,:,:)
|
|---|
| 470 | real(r8), intent(inout) :: avgUT(LBi:,LBj:,:,:)
|
|---|
| 471 | real(r8), intent(inout) :: avgVT(LBi:,LBj:,:,:)
|
|---|
| 472 | real(r8), intent(inout) :: avgTT(LBi:,LBj:,:,:)
|
|---|
| 473 | # endif
|
|---|
| 474 | real(r8), intent(inout) :: avgU2(LBi:,LBj:)
|
|---|
| 475 | real(r8), intent(inout) :: avgV2(LBi:,LBj:)
|
|---|
| 476 | real(r8), intent(inout) :: avgZZ(LBi:,LBj:)
|
|---|
| 477 | # endif
|
|---|
| 478 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 479 | real(r8), intent(inout) :: avgu2Sd(LBi:,LBj:)
|
|---|
| 480 | real(r8), intent(inout) :: avgv2Sd(LBi:,LBj:)
|
|---|
| 481 | real(r8), intent(inout) :: avgu2RS(LBi:,LBj:)
|
|---|
| 482 | real(r8), intent(inout) :: avgv2RS(LBi:,LBj:)
|
|---|
| 483 | real(r8), intent(inout) :: avgSxx2d(LBi:,LBj:)
|
|---|
| 484 | real(r8), intent(inout) :: avgSxy2d(LBi:,LBj:)
|
|---|
| 485 | real(r8), intent(inout) :: avgSyy2d(LBi:,LBj:)
|
|---|
| 486 | # endif
|
|---|
| 487 | real(r8), intent(inout) :: avgu2d(LBi:,LBj:)
|
|---|
| 488 | real(r8), intent(inout) :: avgv2d(LBi:,LBj:)
|
|---|
| 489 | real(r8), intent(inout) :: avgzeta(LBi:,LBj:)
|
|---|
| 490 |
|
|---|
| 491 | # else
|
|---|
| 492 |
|
|---|
| 493 | # ifdef SOLVE3D
|
|---|
| 494 | real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
|
|---|
| 495 | real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
|
|---|
| 496 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 497 | real(r8), intent(in) :: Huon(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 498 | real(r8), intent(in) :: Hvom(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 499 | # endif
|
|---|
| 500 | real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
|
|---|
| 501 | real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
|
|---|
| 502 | real(r8), intent(in) :: W(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 503 | real(r8), intent(in) :: wvel(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 504 | real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
|---|
| 505 | real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 506 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 507 | real(r8), intent(in) :: u_stokes(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 508 | real(r8), intent(in) :: v_stokes(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 509 | real(r8), intent(in) :: rustr3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 510 | real(r8), intent(in) :: rvstr3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 511 | real(r8), intent(in) :: Sxx(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 512 | real(r8), intent(in) :: Sxy(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 513 | real(r8), intent(in) :: Syy(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 514 | real(r8), intent(in) :: Szx(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 515 | real(r8), intent(in) :: Szy(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 516 | # endif
|
|---|
| 517 | # ifdef LMD_SKPP
|
|---|
| 518 | real(r8), intent(in) :: hsbl(LBi:UBi,LBj:UBj)
|
|---|
| 519 | # endif
|
|---|
| 520 | # ifdef LMD_BKPP
|
|---|
| 521 | real(r8), intent(in) :: hbbl(LBi:UBi,LBj:UBj)
|
|---|
| 522 | # endif
|
|---|
| 523 | # ifdef AVERAGES_AKV
|
|---|
| 524 | real(r8), intent(in) :: Akv(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 525 | # endif
|
|---|
| 526 | # if defined AVERAGES_AKT || defined AVERAGES_AKS
|
|---|
| 527 | real(r8), intent(in) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
|
|---|
| 528 | # endif
|
|---|
| 529 | # ifdef AVERAGES_FLUXES
|
|---|
| 530 | real(r8), intent(in) :: stflx(LBi:UBi,LBj:UBj,NT(ng))
|
|---|
| 531 | # ifdef BULK_FLUXES
|
|---|
| 532 | real(r8), intent(in) :: lhflx(LBi:UBi,LBj:UBj)
|
|---|
| 533 | real(r8), intent(in) :: shflx(LBi:UBi,LBj:UBj)
|
|---|
| 534 | real(r8), intent(in) :: lrflx(LBi:UBi,LBj:UBj)
|
|---|
| 535 | # ifdef EMINUSP
|
|---|
| 536 | real(r8), intent(in) :: evap(LBi:UBi,LBj:UBj)
|
|---|
| 537 | real(r8), intent(in) :: rain(LBi:UBi,LBj:UBj)
|
|---|
| 538 | # endif
|
|---|
| 539 | # endif
|
|---|
| 540 | # ifdef SHORTWAVE
|
|---|
| 541 | real(r8), intent(in) :: srflx(LBi:UBi,LBj:UBj)
|
|---|
| 542 | # endif
|
|---|
| 543 | # endif
|
|---|
| 544 | # endif
|
|---|
| 545 | # ifdef AVERAGES_FLUXES
|
|---|
| 546 | real(r8), intent(in) :: sustr(LBi:UBi,LBj:UBj)
|
|---|
| 547 | real(r8), intent(in) :: svstr(LBi:UBi,LBj:UBj)
|
|---|
| 548 | # endif
|
|---|
| 549 | real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,3)
|
|---|
| 550 | real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,3)
|
|---|
| 551 | real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj,3)
|
|---|
| 552 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 553 | real(r8), intent(in) :: ubar_stokes(LBi:UBi,LBj:UBj)
|
|---|
| 554 | real(r8), intent(in) :: vbar_stokes(LBi:UBi,LBj:UBj)
|
|---|
| 555 | real(r8), intent(in) :: rustr2d(LBi:UBi,LBj:UBj)
|
|---|
| 556 | real(r8), intent(in) :: rvstr2d(LBi:UBi,LBj:UBj)
|
|---|
| 557 | real(r8), intent(in) :: Sxx_bar(LBi:UBi,LBj:UBj)
|
|---|
| 558 | real(r8), intent(in) :: Sxy_bar(LBi:UBi,LBj:UBj)
|
|---|
| 559 | real(r8), intent(in) :: Syy_bar(LBi:UBi,LBj:UBj)
|
|---|
| 560 | # endif
|
|---|
| 561 | # ifdef SOLVE3D
|
|---|
| 562 | real(r8), intent(inout) :: avgu3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 563 | real(r8), intent(inout) :: avgv3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 564 | real(r8), intent(inout) :: avgw3d(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 565 | real(r8), intent(inout) :: avgwvel(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 566 | real(r8), intent(inout) :: avgt(LBi:UBi,LBj:UBj,N(ng),NT(ng))
|
|---|
| 567 | real(r8), intent(inout) :: avgrho(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 568 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 569 | real(r8), intent(inout) :: avgu3Sd(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 570 | real(r8), intent(inout) :: avgv3Sd(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 571 | real(r8), intent(inout) :: avgu3RS(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 572 | real(r8), intent(inout) :: avgv3RS(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 573 | real(r8), intent(inout) :: avgSxx3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 574 | real(r8), intent(inout) :: avgSxy3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 575 | real(r8), intent(inout) :: avgSyy3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 576 | real(r8), intent(inout) :: avgSzx3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 577 | real(r8), intent(inout) :: avgSzy3d(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 578 | # endif
|
|---|
| 579 | # ifdef LMD_SKPP
|
|---|
| 580 | real(r8), intent(inout) :: avghsbl(LBi:UBi,LBj:UBj)
|
|---|
| 581 | # endif
|
|---|
| 582 | # ifdef LMD_BKPP
|
|---|
| 583 | real(r8), intent(inout) :: avghbbl(LBi:UBi,LBj:UBj)
|
|---|
| 584 | # endif
|
|---|
| 585 | # ifdef AVERAGES_AKV
|
|---|
| 586 | real(r8), intent(inout) :: avgAKv(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 587 | # endif
|
|---|
| 588 | # ifdef AVERAGES_AKT
|
|---|
| 589 | real(r8), intent(inout) :: avgAKt(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 590 | # endif
|
|---|
| 591 | # ifdef AVERAGES_AKS
|
|---|
| 592 | real(r8), intent(inout) :: avgAKs(LBi:UBi,LBj:UBj,0:N(ng))
|
|---|
| 593 | # endif
|
|---|
| 594 | # ifdef AVERAGES_FLUXES
|
|---|
| 595 | real(r8), intent(inout) :: avgstf(LBi:UBi,LBj:UBj)
|
|---|
| 596 | real(r8), intent(inout) :: avgswf(LBi:UBi,LBj:UBj)
|
|---|
| 597 | # ifdef BULK_FLUXES
|
|---|
| 598 | real(r8), intent(inout) :: avglhf(LBi:UBi,LBj:UBj)
|
|---|
| 599 | real(r8), intent(inout) :: avgshf(LBi:UBi,LBj:UBj)
|
|---|
| 600 | real(r8), intent(inout) :: avglrf(LBi:UBi,LBj:UBj)
|
|---|
| 601 | # ifdef EMINUSP
|
|---|
| 602 | real(r8), intent(inout) :: avgevap(LBi:UBi,LBj:UBj)
|
|---|
| 603 | real(r8), intent(inout) :: avgrain(LBi:UBi,LBj:UBj)
|
|---|
| 604 | # endif
|
|---|
| 605 | # endif
|
|---|
| 606 | # ifdef SHORTWAVE
|
|---|
| 607 | real(r8), intent(inout) :: avgsrf(LBi:UBi,LBj:UBj)
|
|---|
| 608 | # endif
|
|---|
| 609 | # endif
|
|---|
| 610 | # endif
|
|---|
| 611 | # ifdef AVERAGES_FLUXES
|
|---|
| 612 | real(r8), intent(inout) :: avgsus(LBi:UBi,LBj:UBj)
|
|---|
| 613 | real(r8), intent(inout) :: avgsvs(LBi:UBi,LBj:UBj)
|
|---|
| 614 | # endif
|
|---|
| 615 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 616 | # ifdef SOLVE3D
|
|---|
| 617 | real(r8), intent(inout) :: avgHuon(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 618 | real(r8), intent(inout) :: avgHvom(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 619 | real(r8), intent(inout) :: avgUU(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 620 | real(r8), intent(inout) :: avgUV(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 621 | real(r8), intent(inout) :: avgVV(LBi:UBi,LBj:UBj,N(ng))
|
|---|
| 622 | real(r8), intent(inout) :: avgHuonT(LBi:UBi,LBj:UBj,N(ng),NAT)
|
|---|
| 623 | real(r8), intent(inout) :: avgHvomT(LBi:UBi,LBj:UBj,N(ng),NAT)
|
|---|
| 624 | real(r8), intent(inout) :: avgUT(LBi:UBi,LBj:UBj,N(ng),NAT)
|
|---|
| 625 | real(r8), intent(inout) :: avgVT(LBi:UBi,LBj:UBj,N(ng),NAT)
|
|---|
| 626 | real(r8), intent(inout) :: avgTT(LBi:UBi,LBj:UBj,N(ng),NAT)
|
|---|
| 627 | # endif
|
|---|
| 628 | real(r8), intent(inout) :: avgU2(LBi:UBi,LBj:UBj)
|
|---|
| 629 | real(r8), intent(inout) :: avgV2(LBi:UBi,LBj:UBj)
|
|---|
| 630 | real(r8), intent(inout) :: avgZZ(LBi:UBi,LBj:UBj)
|
|---|
| 631 | # endif
|
|---|
| 632 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 633 | real(r8), intent(inout) :: avgu2Sd(LBi:UBi,LBj:UBj)
|
|---|
| 634 | real(r8), intent(inout) :: avgv2Sd(LBi:UBi,LBj:UBj)
|
|---|
| 635 | real(r8), intent(inout) :: avgu2RS(LBi:UBi,LBj:UBj)
|
|---|
| 636 | real(r8), intent(inout) :: avgv2RS(LBi:UBi,LBj:UBj)
|
|---|
| 637 | real(r8), intent(inout) :: avgSxx2d(LBi:UBi,LBj:UBj)
|
|---|
| 638 | real(r8), intent(inout) :: avgSxy2d(LBi:UBi,LBj:UBj)
|
|---|
| 639 | real(r8), intent(inout) :: avgSyy2d(LBi:UBi,LBj:UBj)
|
|---|
| 640 | # endif
|
|---|
| 641 | real(r8), intent(inout) :: avgu2d(LBi:UBi,LBj:UBj)
|
|---|
| 642 | real(r8), intent(inout) :: avgv2d(LBi:UBi,LBj:UBj)
|
|---|
| 643 | real(r8), intent(inout) :: avgzeta(LBi:UBi,LBj:UBj)
|
|---|
| 644 | # endif
|
|---|
| 645 | !
|
|---|
| 646 | ! Local variable declarations.
|
|---|
| 647 | !
|
|---|
| 648 | integer :: IstrR, IendR, JstrR, JendR, IstrU, JstrV
|
|---|
| 649 | integer :: i, itrc, j, k
|
|---|
| 650 |
|
|---|
| 651 | real(r8) :: fac
|
|---|
| 652 |
|
|---|
| 653 | # include "set_bounds.h"
|
|---|
| 654 | !
|
|---|
| 655 | !-----------------------------------------------------------------------
|
|---|
| 656 | ! Initialize time-averaged arrays when appropriate. Notice that
|
|---|
| 657 | ! fields are initilized twice during re-start. However, the time-
|
|---|
| 658 | ! averaged fields are computed correctly.
|
|---|
| 659 | !-----------------------------------------------------------------------
|
|---|
| 660 | !
|
|---|
| 661 | IF (((iic(ng).gt.ntsAVG(ng)).and. &
|
|---|
| 662 | & (MOD(iic(ng)-1,nAVG(ng)).eq.1)).or. &
|
|---|
| 663 | & ((nrrec(ng).gt.0).and.(iic(ng).eq.ntstart(ng)))) THEN
|
|---|
| 664 | !
|
|---|
| 665 | ! Initialize 2D fields.
|
|---|
| 666 | !
|
|---|
| 667 | DO j=JstrR,JendR
|
|---|
| 668 | DO i=IstrR,IendR
|
|---|
| 669 | avgzeta(i,j)=zeta(i,j,Kout)
|
|---|
| 670 | avgu2d (i,j)=ubar(i,j,Kout)
|
|---|
| 671 | avgv2d (i,j)=vbar(i,j,Kout)
|
|---|
| 672 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 673 | avgu2Sd(i,j)=ubar_stokes(i,j)
|
|---|
| 674 | avgv2Sd(i,j)=vbar_stokes(i,j)
|
|---|
| 675 | avgu2RS(i,j)=rustr2d(i,j)
|
|---|
| 676 | avgv2RS(i,j)=rvstr2d(i,j)
|
|---|
| 677 | avgSxx2d(i,j)=Sxx_bar(i,j)
|
|---|
| 678 | avgSxy2d(i,j)=Sxy_bar(i,j)
|
|---|
| 679 | avgSyy2d(i,j)=Syy_bar(i,j)
|
|---|
| 680 | # endif
|
|---|
| 681 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 682 | avgZZ(i,j)=zeta(i,j,Kout)*zeta(i,j,Kout)
|
|---|
| 683 | avgU2(i,j)=ubar(i,j,Kout)*ubar(i,j,Kout)
|
|---|
| 684 | avgV2(i,j)=vbar(i,j,Kout)*vbar(i,j,Kout)
|
|---|
| 685 | # endif
|
|---|
| 686 | # ifdef SOLVE3D
|
|---|
| 687 | # ifdef LMD_SKPP
|
|---|
| 688 | avghsbl(i,j)=hsbl(i,j)
|
|---|
| 689 | # endif
|
|---|
| 690 | # ifdef LMD_BKPP
|
|---|
| 691 | avghbbl(i,j)=hbbl(i,j)
|
|---|
| 692 | # endif
|
|---|
| 693 | # ifdef AVERAGES_FLUXES
|
|---|
| 694 | avgstf(i,j)=stflx(i,j,itemp)
|
|---|
| 695 | avgswf(i,j)=stflx(i,j,isalt)
|
|---|
| 696 | # ifdef BULK_FLUXES
|
|---|
| 697 | avglhf(i,j)=lhflx(i,j)
|
|---|
| 698 | avgshf(i,j)=shflx(i,j)
|
|---|
| 699 | avglrf(i,j)=lrflx(i,j)
|
|---|
| 700 | # ifdef EMINUSP
|
|---|
| 701 | avgevap(i,j)=evap(i,j)
|
|---|
| 702 | avgrain(i,j)=rain(i,j)
|
|---|
| 703 | # endif
|
|---|
| 704 | # endif
|
|---|
| 705 | # ifdef SHORTWAVE
|
|---|
| 706 | avgsrf(i,j)=srflx(i,j)
|
|---|
| 707 | # endif
|
|---|
| 708 | # endif
|
|---|
| 709 | # endif
|
|---|
| 710 | # ifdef AVERAGES_FLUXES
|
|---|
| 711 | avgsus(i,j)=sustr(i,j)
|
|---|
| 712 | avgsvs(i,j)=svstr(i,j)
|
|---|
| 713 | # endif
|
|---|
| 714 | END DO
|
|---|
| 715 | END DO
|
|---|
| 716 |
|
|---|
| 717 | # ifdef SOLVE3D
|
|---|
| 718 | !
|
|---|
| 719 | ! Initialize fields associated with 3D horizontal momentum.
|
|---|
| 720 | !
|
|---|
| 721 | DO k=1,N(ng)
|
|---|
| 722 | DO j=JstrR,JendR
|
|---|
| 723 | DO i=IstrR,IendR
|
|---|
| 724 | avgu3d(i,j,k)=u(i,j,k,Nout)
|
|---|
| 725 | avgv3d(i,j,k)=v(i,j,k,Nout)
|
|---|
| 726 | avgrho(i,j,k)=rho(i,j,k)
|
|---|
| 727 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 728 | avgu3Sd(i,j,k)=u_stokes(i,j,k)
|
|---|
| 729 | avgv3Sd(i,j,k)=v_stokes(i,j,k)
|
|---|
| 730 | avgu3RS(i,j,k)=rustr3d(i,j,k)
|
|---|
| 731 | avgv3RS(i,j,k)=rvstr3d(i,j,k)
|
|---|
| 732 | avgSxx3d(i,j,k)=Sxx(i,j,k)
|
|---|
| 733 | avgSxy3d(i,j,k)=Sxy(i,j,k)
|
|---|
| 734 | avgSyy3d(i,j,k)=Syy(i,j,k)
|
|---|
| 735 | avgSzx3d(i,j,k)=Szx(i,j,k)
|
|---|
| 736 | avgSzy3d(i,j,k)=Szy(i,j,k)
|
|---|
| 737 | # endif
|
|---|
| 738 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 739 | avgHuon(i,j,k)=Huon(i,j,k)
|
|---|
| 740 | avgHvom(i,j,k)=Hvom(i,j,k)
|
|---|
| 741 | avgUU(i,j,k)=u(i,j,k,Nout)*u(i,j,k,Nout)
|
|---|
| 742 | avgVV(i,j,k)=v(i,j,k,Nout)*v(i,j,k,Nout)
|
|---|
| 743 | # endif
|
|---|
| 744 | END DO
|
|---|
| 745 | END DO
|
|---|
| 746 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 747 | DO j=Jstr,Jend
|
|---|
| 748 | DO i=Istr,Iend
|
|---|
| 749 | avgUV(i,j,k)=0.25_r8*(u(i ,j ,k,Nout)+ &
|
|---|
| 750 | & u(i+1,j ,k,Nout))* &
|
|---|
| 751 | & (v(i ,j ,k,Nout)+ &
|
|---|
| 752 | & v(i ,j+1,k,Nout))
|
|---|
| 753 | END DO
|
|---|
| 754 | END DO
|
|---|
| 755 | # endif
|
|---|
| 756 | END DO
|
|---|
| 757 | !
|
|---|
| 758 | ! Initialized fields associated with tracers.
|
|---|
| 759 | !
|
|---|
| 760 | DO itrc=1,NT(ng)
|
|---|
| 761 | DO k=1,N(ng)
|
|---|
| 762 | DO j=JstrR,JendR
|
|---|
| 763 | DO i=IstrR,IendR
|
|---|
| 764 | avgt(i,j,k,itrc)=t(i,j,k,Nout,itrc)
|
|---|
| 765 | END DO
|
|---|
| 766 | END DO
|
|---|
| 767 | END DO
|
|---|
| 768 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 769 | IF (itrc.le.NAT) THEN
|
|---|
| 770 | DO k=1,N(ng)
|
|---|
| 771 | DO j=JstrR,JendR
|
|---|
| 772 | DO i=IstrR,IendR
|
|---|
| 773 | avgTT(i,j,k,itrc)=t(i,j,k,Nout,itrc)* &
|
|---|
| 774 | & t(i,j,k,Nout,itrc)
|
|---|
| 775 | END DO
|
|---|
| 776 | DO i=Istr,Iend
|
|---|
| 777 | avgUT(i,j,k,itrc)=0.5_r8*u(i,j,k,Nout)* &
|
|---|
| 778 | & (t(i-1,j,k,Nout,itrc)+ &
|
|---|
| 779 | & t(i ,j,k,Nout,itrc))
|
|---|
| 780 | avgHuonT(i,j,k,itrc)=0.5_r8*Huon(i,j,k)* &
|
|---|
| 781 | & (t(i-1,j,k,Nout,itrc)+ &
|
|---|
| 782 | & t(i ,j,k,Nout,itrc))
|
|---|
| 783 | END DO
|
|---|
| 784 | END DO
|
|---|
| 785 | DO j=Jstr,Jend
|
|---|
| 786 | DO i=IstrR,IendR
|
|---|
| 787 | avgVT(i,j,k,itrc)=0.5_r8*v(i,j,k,Nout)* &
|
|---|
| 788 | & (t(i,j-1,k,Nout,itrc)+ &
|
|---|
| 789 | & t(i,j ,k,Nout,itrc))
|
|---|
| 790 | avgHvomT(i,j,k,itrc)=0.5_r8*Hvom(i,j,k)* &
|
|---|
| 791 | & (t(i,j-1,k,Nout,itrc)+ &
|
|---|
| 792 | & t(i,j ,k,Nout,itrc))
|
|---|
| 793 | END DO
|
|---|
| 794 | END DO
|
|---|
| 795 | END DO
|
|---|
| 796 | END IF
|
|---|
| 797 | # endif
|
|---|
| 798 | END DO
|
|---|
| 799 | !
|
|---|
| 800 | ! Initialize fields at W-points.
|
|---|
| 801 | !
|
|---|
| 802 | DO k=0,N(ng)
|
|---|
| 803 | DO j=JstrR,JendR
|
|---|
| 804 | DO i=IstrR,IendR
|
|---|
| 805 | avgw3d(i,j,k)=W(i,j,k)*pm(i,j)*pn(i,j)
|
|---|
| 806 | avgwvel(i,j,k)=wvel(i,j,k)
|
|---|
| 807 | # ifdef AVERAGES_AKV
|
|---|
| 808 | avgAKv(i,j,k)=Akv(i,j,k)
|
|---|
| 809 | # endif
|
|---|
| 810 | # ifdef AVERAGES_AKT
|
|---|
| 811 | avgAKt(i,j,k)=Akt(i,j,k,itemp)
|
|---|
| 812 | # endif
|
|---|
| 813 | # ifdef AVERAGES_AKS
|
|---|
| 814 | avgAKs(i,j,k)=Akt(i,j,k,isalt)
|
|---|
| 815 | # endif
|
|---|
| 816 | END DO
|
|---|
| 817 | END DO
|
|---|
| 818 | END DO
|
|---|
| 819 | # endif
|
|---|
| 820 | !
|
|---|
| 821 | !-----------------------------------------------------------------------
|
|---|
| 822 | ! Accumulate time-averaged fields.
|
|---|
| 823 | !-----------------------------------------------------------------------
|
|---|
| 824 | !
|
|---|
| 825 | ELSE IF (iic(ng).gt.ntsAVG(ng)) THEN
|
|---|
| 826 | !
|
|---|
| 827 | ! Accumulate 2D fields.
|
|---|
| 828 | !
|
|---|
| 829 | DO j=JstrR,JendR
|
|---|
| 830 | DO i=IstrR,IendR
|
|---|
| 831 | avgzeta(i,j)=avgzeta(i,j)+zeta(i,j,Kout)
|
|---|
| 832 | avgu2d (i,j)=avgu2d (i,j)+ubar(i,j,Kout)
|
|---|
| 833 | avgv2d (i,j)=avgv2d (i,j)+vbar(i,j,Kout)
|
|---|
| 834 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 835 | avgu2Sd(i,j)=avgu2Sd(i,j)+ubar_stokes(i,j)
|
|---|
| 836 | avgv2Sd(i,j)=avgv2Sd(i,j)+vbar_stokes(i,j)
|
|---|
| 837 | avgu2RS(i,j)=avgu2RS(i,j)+rustr2d(i,j)
|
|---|
| 838 | avgv2RS(i,j)=avgv2RS(i,j)+rvstr2d(i,j)
|
|---|
| 839 | avgSxx2d(i,j)=avgSxx2d(i,j)+Sxx_bar(i,j)
|
|---|
| 840 | avgSxy2d(i,j)=avgSxy2d(i,j)+Sxy_bar(i,j)
|
|---|
| 841 | avgSyy2d(i,j)=avgSyy2d(i,j)+Syy_bar(i,j)
|
|---|
| 842 | # endif
|
|---|
| 843 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 844 | avgZZ(i,j)=avgZZ(i,j)+zeta(i,j,Kout)*zeta(i,j,Kout)
|
|---|
| 845 | avgU2(i,j)=avgU2(i,j)+ubar(i,j,Kout)*ubar(i,j,Kout)
|
|---|
| 846 | avgV2(i,j)=avgU2(i,j)+vbar(i,j,Kout)*vbar(i,j,Kout)
|
|---|
| 847 | # endif
|
|---|
| 848 | # ifdef SOLVE3D
|
|---|
| 849 | # ifdef LMD_SKPP
|
|---|
| 850 | avghsbl(i,j)=avghsbl(i,j)+hsbl(i,j)
|
|---|
| 851 | # endif
|
|---|
| 852 | # ifdef LMD_BKPP
|
|---|
| 853 | avghbbl(i,j)=avghbbl(i,j)+hbbl(i,j)
|
|---|
| 854 | # endif
|
|---|
| 855 | # ifdef AVERAGES_FLUXES
|
|---|
| 856 | avgstf(i,j)=avgstf(i,j)+stflx(i,j,itemp)
|
|---|
| 857 | avgswf(i,j)=avgswf(i,j)+stflx(i,j,isalt)
|
|---|
| 858 | # ifdef BULK_FLUXES
|
|---|
| 859 | avglhf(i,j)=avglhf(i,j)+lhflx(i,j)
|
|---|
| 860 | avgshf(i,j)=avgshf(i,j)+shflx(i,j)
|
|---|
| 861 | avglrf(i,j)=avglrf(i,j)+lrflx(i,j)
|
|---|
| 862 | # ifdef EMINUSP
|
|---|
| 863 | avgevap(i,j)=avgevap(i,j)+evap(i,j)
|
|---|
| 864 | avgrain(i,j)=avgrain(i,j)+rain(i,j)
|
|---|
| 865 | # endif
|
|---|
| 866 | # endif
|
|---|
| 867 | # ifdef SHORTWAVE
|
|---|
| 868 | avgsrf(i,j)=avgsrf(i,j)+srflx(i,j)
|
|---|
| 869 | # endif
|
|---|
| 870 | # endif
|
|---|
| 871 | # endif
|
|---|
| 872 | # ifdef AVERAGES_FLUXES
|
|---|
| 873 | avgsus(i,j)=avgsus(i,j)+sustr(i,j)
|
|---|
| 874 | avgsvs(i,j)=avgsvs(i,j)+svstr(i,j)
|
|---|
| 875 | # endif
|
|---|
| 876 | END DO
|
|---|
| 877 | END DO
|
|---|
| 878 |
|
|---|
| 879 | # ifdef SOLVE3D
|
|---|
| 880 | !
|
|---|
| 881 | ! Accumulate fields associated with 3D horizontal momentum.
|
|---|
| 882 | !
|
|---|
| 883 | DO k=1,N(ng)
|
|---|
| 884 | DO j=JstrR,JendR
|
|---|
| 885 | DO i=IstrR,IendR
|
|---|
| 886 | avgu3d(i,j,k)=avgu3d(i,j,k)+u(i,j,k,Nout)
|
|---|
| 887 | avgv3d(i,j,k)=avgv3d(i,j,k)+v(i,j,k,Nout)
|
|---|
| 888 | avgrho(i,j,k)=avgrho(i,j,k)+rho(i,j,k)
|
|---|
| 889 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 890 | avgu3Sd(i,j,k)=avgu3Sd(i,j,k)+u_stokes(i,j,k)
|
|---|
| 891 | avgv3Sd(i,j,k)=avgv3Sd(i,j,k)+v_stokes(i,j,k)
|
|---|
| 892 | avgu3RS(i,j,k)=avgu3RS(i,j,k)+rustr3d(i,j,k)
|
|---|
| 893 | avgv3RS(i,j,k)=avgv3RS(i,j,k)+rvstr3d(i,j,k)
|
|---|
| 894 | avgSxx3d(i,j,k)=avgSxx3d(i,j,k)+Sxx(i,j,k)
|
|---|
| 895 | avgSxy3d(i,j,k)=avgSxy3d(i,j,k)+Sxy(i,j,k)
|
|---|
| 896 | avgSyy3d(i,j,k)=avgSyy3d(i,j,k)+Syy(i,j,k)
|
|---|
| 897 | avgSzx3d(i,j,k)=avgSzx3d(i,j,k)+Szx(i,j,k)
|
|---|
| 898 | avgSzy3d(i,j,k)=avgSzy3d(i,j,k)+Szy(i,j,k)
|
|---|
| 899 | # endif
|
|---|
| 900 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 901 | avgHuon(i,j,k)=avgHuon(i,j,k)+Huon(i,j,k)
|
|---|
| 902 | avgHvom(i,j,k)=avgHvom(i,j,k)+Hvom(i,j,k)
|
|---|
| 903 | avgUU(i,j,k)=avgUU(i,j,k)+u(i,j,k,Nout)*u(i,j,k,Nout)
|
|---|
| 904 | avgVV(i,j,k)=avgVV(i,j,k)+v(i,j,k,Nout)*v(i,j,k,Nout)
|
|---|
| 905 | # endif
|
|---|
| 906 | END DO
|
|---|
| 907 | END DO
|
|---|
| 908 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 909 | DO j=Jstr,Jend
|
|---|
| 910 | DO i=Istr,Iend
|
|---|
| 911 | avgUV(i,j,k)=avgUV(i,j,k)+ &
|
|---|
| 912 | & 0.25_r8*(u(i ,j ,k,Nout)+ &
|
|---|
| 913 | & u(i+1,j ,k,Nout))* &
|
|---|
| 914 | & (v(i ,j ,k,Nout)+ &
|
|---|
| 915 | & v(i ,j+1,k,Nout))
|
|---|
| 916 | END DO
|
|---|
| 917 | END DO
|
|---|
| 918 | # endif
|
|---|
| 919 | END DO
|
|---|
| 920 | !
|
|---|
| 921 | ! Accumulate fields associated with tracers.
|
|---|
| 922 | !
|
|---|
| 923 | DO itrc=1,NT(ng)
|
|---|
| 924 | DO k=1,N(ng)
|
|---|
| 925 | DO j=JstrR,JendR
|
|---|
| 926 | DO i=IstrR,IendR
|
|---|
| 927 | avgt(i,j,k,itrc)=avgt(i,j,k,itrc)+t(i,j,k,Nout,itrc)
|
|---|
| 928 | END DO
|
|---|
| 929 | END DO
|
|---|
| 930 | END DO
|
|---|
| 931 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 932 | IF (itrc.le.NAT) THEN
|
|---|
| 933 | DO k=1,N(ng)
|
|---|
| 934 | DO j=JstrR,JendR
|
|---|
| 935 | DO i=IstrR,IendR
|
|---|
| 936 | avgTT(i,j,k,itrc)=avgTT(i,j,k,itrc)+ &
|
|---|
| 937 | & t(i,j,k,Nout,itrc)* &
|
|---|
| 938 | & t(i,j,k,Nout,itrc)
|
|---|
| 939 | END DO
|
|---|
| 940 | DO i=Istr,Iend
|
|---|
| 941 | avgUT(i,j,k,itrc)=avgUT(i,j,k,itrc)+ &
|
|---|
| 942 | & 0.5_r8*u(i,j,k,Nout)* &
|
|---|
| 943 | & (t(i-1,j,k,Nout,itrc)+ &
|
|---|
| 944 | & t(i ,j,k,Nout,itrc))
|
|---|
| 945 | avgHuonT(i,j,k,itrc)=avgHuonT(i,j,k,itrc)+ &
|
|---|
| 946 | & 0.5_r8*Huon(i,j,k)* &
|
|---|
| 947 | & (t(i-1,j,k,Nout,itrc)+ &
|
|---|
| 948 | & t(i ,j,k,Nout,itrc))
|
|---|
| 949 | END DO
|
|---|
| 950 | END DO
|
|---|
| 951 | DO j=Jstr,Jend
|
|---|
| 952 | DO i=IstrR,IendR
|
|---|
| 953 | avgVT(i,j,k,itrc)=avgVT(i,j,k,itrc)+ &
|
|---|
| 954 | & 0.5_r8*v(i,j,k,Nout)* &
|
|---|
| 955 | & (t(i,j-1,k,Nout,itrc)+ &
|
|---|
| 956 | & t(i,j ,k,Nout,itrc))
|
|---|
| 957 | avgHvomT(i,j,k,itrc)=avgVT(i,j,k,itrc)+ &
|
|---|
| 958 | & 0.5_r8*Hvom(i,j,k)* &
|
|---|
| 959 | & (t(i,j-1,k,Nout,itrc)+ &
|
|---|
| 960 | & t(i,j ,k,Nout,itrc))
|
|---|
| 961 | END DO
|
|---|
| 962 | END DO
|
|---|
| 963 | END DO
|
|---|
| 964 | END IF
|
|---|
| 965 | # endif
|
|---|
| 966 | END DO
|
|---|
| 967 | !
|
|---|
| 968 | ! Accumulate fields at W-points.
|
|---|
| 969 | !
|
|---|
| 970 | DO k=0,N(ng)
|
|---|
| 971 | DO j=JstrR,JendR
|
|---|
| 972 | DO i=IstrR,IendR
|
|---|
| 973 | avgw3d(i,j,k)=avgw3d(i,j,k)+W(i,j,k)*pm(i,j)*pn(i,j)
|
|---|
| 974 | avgwvel(i,j,k)=avgwvel(i,j,k)+wvel(i,j,k)
|
|---|
| 975 | # ifdef AVERAGES_AKV
|
|---|
| 976 | avgAKv(i,j,k)=avgAKv(i,j,k)+Akv(i,j,k)
|
|---|
| 977 | # endif
|
|---|
| 978 | # ifdef AVERAGES_AKT
|
|---|
| 979 | avgAKt(i,j,k)=avgAKt(i,j,k)+Akt(i,j,k,itemp)
|
|---|
| 980 | # endif
|
|---|
| 981 | # ifdef AVERAGES_AKS
|
|---|
| 982 | avgAKs(i,j,k)=avgAKs(i,j,k)+Akt(i,j,k,isalt)
|
|---|
| 983 | # endif
|
|---|
| 984 | END DO
|
|---|
| 985 | END DO
|
|---|
| 986 | END DO
|
|---|
| 987 | # endif
|
|---|
| 988 | END IF
|
|---|
| 989 | !
|
|---|
| 990 | !-----------------------------------------------------------------------
|
|---|
| 991 | ! Convert accumulated sums into time-averages, if appropriate.
|
|---|
| 992 | !-----------------------------------------------------------------------
|
|---|
| 993 | !
|
|---|
| 994 | IF ((iic(ng).gt.ntsAVG(ng)).and. &
|
|---|
| 995 | & (MOD(iic(ng)-1,nAVG(ng)).eq.0).and. &
|
|---|
| 996 | & ((iic(ng).ne.ntstart(ng)).or.(nrrec(ng).eq.0))) THEN
|
|---|
| 997 | fac=1.0_r8/REAL(nAVG(ng),r8)
|
|---|
| 998 | IF (SOUTH_WEST_TEST) THEN
|
|---|
| 999 | AVGtime(ng)=AVGtime(ng)+REAL(nAVG(ng),r8)*dt(ng)
|
|---|
| 1000 | END IF
|
|---|
| 1001 | !
|
|---|
| 1002 | ! Process 2D fields.
|
|---|
| 1003 | !
|
|---|
| 1004 | DO j=JstrR,JendR
|
|---|
| 1005 | DO i=IstrR,IendR
|
|---|
| 1006 | avgzeta(i,j)=fac*avgzeta(i,j)
|
|---|
| 1007 | avgu2d (i,j)=fac*avgu2d (i,j)
|
|---|
| 1008 | avgv2d (i,j)=fac*avgv2d (i,j)
|
|---|
| 1009 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 1010 | avgu2Sd(i,j)=fac*avgu2Sd(i,j)
|
|---|
| 1011 | avgv2Sd(i,j)=fac*avgv2Sd(i,j)
|
|---|
| 1012 | avgu2RS(i,j)=fac*avgu2RS(i,j)
|
|---|
| 1013 | avgv2RS(i,j)=fac*avgv2RS(i,j)
|
|---|
| 1014 | avgSxx2d(i,j)=fac*avgSxx2d(i,j)
|
|---|
| 1015 | avgSxy2d(i,j)=fac*avgSxy2d(i,j)
|
|---|
| 1016 | avgSyy2d(i,j)=fac*avgSyy2d(i,j)
|
|---|
| 1017 | # endif
|
|---|
| 1018 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 1019 | avgZZ(i,j)=fac*avgZZ(i,j)
|
|---|
| 1020 | avgU2(i,j)=fac*avgU2(i,j)
|
|---|
| 1021 | avgV2(i,j)=fac*avgU2(i,j)
|
|---|
| 1022 | # endif
|
|---|
| 1023 | # ifdef SOLVE3D
|
|---|
| 1024 | # ifdef LMD_SKPP
|
|---|
| 1025 | avghsbl(i,j)=fac*avghsbl(i,j)
|
|---|
| 1026 | # endif
|
|---|
| 1027 | # ifdef LMD_BKPP
|
|---|
| 1028 | avghbbl(i,j)=fac*avghbbl(i,j)
|
|---|
| 1029 | # endif
|
|---|
| 1030 | # ifdef AVERAGES_FLUXES
|
|---|
| 1031 | avgstf(i,j)=fac*avgstf(i,j)
|
|---|
| 1032 | avgswf(i,j)=fac*avgswf(i,j)
|
|---|
| 1033 | # ifdef BULK_FLUXES
|
|---|
| 1034 | avglhf(i,j)=fac*avglhf(i,j)
|
|---|
| 1035 | avgshf(i,j)=fac*avgshf(i,j)
|
|---|
| 1036 | avglrf(i,j)=fac*avglrf(i,j)
|
|---|
| 1037 | # ifdef EMINUSP
|
|---|
| 1038 | avgevap(i,j)=fac*avgevap(i,j)
|
|---|
| 1039 | avgrain(i,j)=fac*avgrain(i,j)
|
|---|
| 1040 | # endif
|
|---|
| 1041 | # endif
|
|---|
| 1042 | # ifdef SHORTWAVE
|
|---|
| 1043 | avgsrf(i,j)=fac*avgsrf(i,j)
|
|---|
| 1044 | # endif
|
|---|
| 1045 | # endif
|
|---|
| 1046 | # endif
|
|---|
| 1047 | # ifdef AVERAGES_FLUXES
|
|---|
| 1048 | avgsus(i,j)=fac*avgsus(i,j)
|
|---|
| 1049 | avgsvs(i,j)=fac*avgsvs(i,j)
|
|---|
| 1050 | # endif
|
|---|
| 1051 | END DO
|
|---|
| 1052 | END DO
|
|---|
| 1053 |
|
|---|
| 1054 | # ifdef SOLVE3D
|
|---|
| 1055 | !
|
|---|
| 1056 | ! Process fields associated with 3D horizontal momentum.
|
|---|
| 1057 | !
|
|---|
| 1058 | DO k=1,N(ng)
|
|---|
| 1059 | DO j=JstrR,JendR
|
|---|
| 1060 | DO i=IstrR,IendR
|
|---|
| 1061 | avgu3d(i,j,k)=fac*avgu3d(i,j,k)
|
|---|
| 1062 | avgv3d(i,j,k)=fac*avgv3d(i,j,k)
|
|---|
| 1063 | avgrho(i,j,k)=fac*avgrho(i,j,k)
|
|---|
| 1064 | # ifdef AVERAGES_NEARSHORE
|
|---|
| 1065 | avgu3Sd(i,j,k)=fac*avgu3Sd(i,j,k)
|
|---|
| 1066 | avgv3Sd(i,j,k)=fac*avgv3Sd(i,j,k)
|
|---|
| 1067 | avgu3RS(i,j,k)=fac*avgu3RS(i,j,k)
|
|---|
| 1068 | avgv3RS(i,j,k)=fac*avgv3RS(i,j,k)
|
|---|
| 1069 | avgSxx3d(i,j,k)=fac*avgSxx3d(i,j,k)
|
|---|
| 1070 | avgSxy3d(i,j,k)=fac*avgSxy3d(i,j,k)
|
|---|
| 1071 | avgSyy3d(i,j,k)=fac*avgSyy3d(i,j,k)
|
|---|
| 1072 | avgSzx3d(i,j,k)=fac*avgSzx3d(i,j,k)
|
|---|
| 1073 | avgSzy3d(i,j,k)=fac*avgSzy3d(i,j,k)
|
|---|
| 1074 | # endif
|
|---|
| 1075 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 1076 | avgHuon(i,j,k)=fac*avgHuon(i,j,k)
|
|---|
| 1077 | avgHvom(i,j,k)=fac*avgHvom(i,j,k)
|
|---|
| 1078 | avgUU(i,j,k)=fac*avgUU(i,j,k)
|
|---|
| 1079 | avgVV(i,j,k)=fac*avgVV(i,j,k)
|
|---|
| 1080 | # endif
|
|---|
| 1081 | END DO
|
|---|
| 1082 | END DO
|
|---|
| 1083 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 1084 | DO j=Jstr,Jend
|
|---|
| 1085 | DO i=Istr,Iend
|
|---|
| 1086 | avgUV(i,j,k)=fac*avgUV(i,j,k)
|
|---|
| 1087 | END DO
|
|---|
| 1088 | END DO
|
|---|
| 1089 | # endif
|
|---|
| 1090 | END DO
|
|---|
| 1091 | !
|
|---|
| 1092 | ! Process fields associated with tracers.
|
|---|
| 1093 | !
|
|---|
| 1094 | DO itrc=1,NT(ng)
|
|---|
| 1095 | DO k=1,N(ng)
|
|---|
| 1096 | DO j=JstrR,JendR
|
|---|
| 1097 | DO i=IstrR,IendR
|
|---|
| 1098 | avgt(i,j,k,itrc)=fac*avgt(i,j,k,itrc)
|
|---|
| 1099 | END DO
|
|---|
| 1100 | END DO
|
|---|
| 1101 | END DO
|
|---|
| 1102 | # ifdef AVERAGES_QUADRATIC
|
|---|
| 1103 | IF (itrc.le.NAT) THEN
|
|---|
| 1104 | DO k=1,N(ng)
|
|---|
| 1105 | DO j=JstrR,JendR
|
|---|
| 1106 | DO i=IstrR,IendR
|
|---|
| 1107 | avgTT(i,j,k,itrc)=fac*avgTT(i,j,k,itrc)
|
|---|
| 1108 | END DO
|
|---|
| 1109 | DO i=Istr,Iend
|
|---|
| 1110 | avgUT(i,j,k,itrc)=fac*avgUT(i,j,k,itrc)
|
|---|
| 1111 | avgHuonT(i,j,k,itrc)=fac*avgHuonT(i,j,k,itrc)
|
|---|
| 1112 | END DO
|
|---|
| 1113 | END DO
|
|---|
| 1114 | DO j=Jstr,Jend
|
|---|
| 1115 | DO i=IstrR,IendR
|
|---|
| 1116 | avgVT(i,j,k,itrc)=fac*avgVT(i,j,k,itrc)
|
|---|
| 1117 | avgHvomT(i,j,k,itrc)=fac*avgHvomT(i,j,k,itrc)
|
|---|
| 1118 | END DO
|
|---|
| 1119 | END DO
|
|---|
| 1120 | END DO
|
|---|
| 1121 | END IF
|
|---|
| 1122 | # endif
|
|---|
| 1123 | END DO
|
|---|
| 1124 | !
|
|---|
| 1125 | ! Process fields at W-points.
|
|---|
| 1126 | !
|
|---|
| 1127 | DO k=0,N(ng)
|
|---|
| 1128 | DO j=JstrR,JendR
|
|---|
| 1129 | DO i=IstrR,IendR
|
|---|
| 1130 | avgw3d(i,j,k)=fac*avgw3d(i,j,k)
|
|---|
| 1131 | avgwvel(i,j,k)=fac*avgwvel(i,j,k)
|
|---|
| 1132 | # ifdef AVERAGES_AKV
|
|---|
| 1133 | avgAKv(i,j,k)=fac*avgAKv(i,j,k)
|
|---|
| 1134 | # endif
|
|---|
| 1135 | # ifdef AVERAGES_AKT
|
|---|
| 1136 | avgAKt(i,j,k)=fac*avgAKt(i,j,k)
|
|---|
| 1137 | # endif
|
|---|
| 1138 | # ifdef AVERAGES_AKS
|
|---|
| 1139 | avgAKs(i,j,k)=fac*avgAKs(i,j,k)
|
|---|
| 1140 | # endif
|
|---|
| 1141 | END DO
|
|---|
| 1142 | END DO
|
|---|
| 1143 | END DO
|
|---|
| 1144 | # endif
|
|---|
| 1145 | END IF
|
|---|
| 1146 |
|
|---|
| 1147 | RETURN
|
|---|
| 1148 | END SUBROUTINE set_avg_tile
|
|---|
| 1149 | #endif
|
|---|
| 1150 | END MODULE set_avg_mod
|
|---|