/[MITgcm]/MITgcm/pkg/opps/opps_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/opps/opps_calc.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.14 - (show annotations) (download)
Thu Mar 10 20:57:11 2016 UTC (8 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65v, checkpoint65w, checkpoint65u, HEAD
Changes since 1.13: +24 -10 lines
- add options (parameter selectP_inEOS_Zc) to select which pressure to use
  in EOS for height coordinate.

1 C $Header: /u/gcmpack/MITgcm/pkg/opps/opps_calc.F,v 1.13 2016/02/15 18:01:59 jmc Exp $
2 C $Name: $
3
4 #include "OPPS_OPTIONS.h"
5 #undef OPPS_ORGCODE
6
7 C-- File opps_calc.F:
8 C-- Contents
9 C-- o OPPS_CALC
10 C-- o STATE1
11 C-- o NLOPPS
12
13 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
14 CBOP
15 C !ROUTINE: OPPS_CALC
16
17 C !INTERFACE: ======================================================
18 SUBROUTINE OPPS_CALC(
19 U tracerEnv,
20 I wVel,
21 I kMax, nTracer, nTracerInuse,
22 I I, J, bi, bj, myTime, myIter, myThid )
23
24 C !DESCRIPTION: \bv
25 C *=====================================================================*
26 C | SUBROUTINE OPPS_CALC |
27 C | o Compute all OPPS fields defined in OPPS.h |
28 C *=====================================================================*
29 C | This subroutine is based on the routine 3dconvection.F |
30 C | by E. Skyllingstad (?) |
31 C | plenty of modifications to make it work: |
32 C | - removed many unused parameters and variables |
33 C | - turned everything (back) into 1D code |
34 C | - pass variables, that are originally in common blocks: |
35 C | maxDepth |
36 C | - pass vertical velocity, set in OPPS_INTERFACE |
37 C | - do not use convadj for now (whatever that is) |
38 C | - changed two .LT. 0 to .LE. 0 statements (because of possible |
39 C | division) |
40 C | - replaced statement function state1 by call to a real function |
41 C | - removed range check, actually moved it up to OPPS_INTERFACE |
42 C | - avoid division by zero: if (Wd.EQ.0) dt = ...1/Wd |
43 C | - cleaned-up debugging |
44 C | - replaced local dz and GridThickness by global drF |
45 C | - replaced 1/dz by 1*recip_drF |
46 C | - replaced 9.81 with gravity (=9.81) |
47 C | - added a lot of comments that relate code to equation in paper |
48 C | (Paluszkiewicz+Romea, 1997, Dynamics of Atmospheres and Oceans, |
49 C | 26, pp. 95-130) |
50 C | - included passive tracer support. This is the main change and may |
51 C | not improve the readability of the code because of the joint |
52 C | treatment of active (theta, salt) and passive tracers. The array |
53 C | tracerEnv(Nr,2+PTRACERS_num) contains |
54 C | theta = tracerEnv(:,1), |
55 C | salt = tracerEnv(:,2), and |
56 C | ptracers = tracerEnv(:,3:PTRACERS_num+2). |
57 C | All related array names have been changed accordingly, so that |
58 C | instead of Sd(Nr) and Td(Nr) (plume salinity and temperature), we |
59 C | have Pd(Nr,nTracer) (tracer in plume), with Sd(:) = Pd(:,2), |
60 C | Td(:) = Pd(:,1), etc. |
61 C | o TODO: |
62 C | clean up the logic of the vertical loops and get rid off the |
63 C | GOTO statements |
64 C *=====================================================================*
65 C \ev
66
67 IMPLICIT NONE
68
69 C !USES: ============================================================
70 #include "SIZE.h"
71 #include "EEPARAMS.h"
72 #include "PARAMS.h"
73 #include "OPPS.h"
74 #include "GRID.h"
75
76 C !INPUT PARAMETERS: ===================================================
77 C Routine arguments
78 C bi, bj :: array indices on which to apply calculations
79 C myTime :: Current time in simulation
80
81 INTEGER I, J, bi, bj, KMax, nTracer, nTracerInUse
82 INTEGER myThid, myIter
83 _RL myTime
84 _RL tracerEnv(Nr,nTracer),wVel(Nr)
85
86 #ifdef ALLOW_OPPS
87 C !FUNCTIONS: ==========================================================
88 c EXTERNAL DIFFERENT_MULTIPLE
89 c LOGICAL DIFFERENT_MULTIPLE
90 _RL STATE1
91
92 C !LOCAL VARIABLES: ====================================================
93 C Local constants
94 C msgBuf :: Informational/error message buffer
95 CHARACTER*(MAX_LEN_MBUF) msgBuf
96 INTEGER K, K2, K2m1, K2p1, ktr
97 INTEGER ntime,nn,kmx,ic
98 INTEGER maxDepth
99
100 _RL wsqr,oldflux,newflux,entrainrate
101 _RL pmix
102 _RL D1, D2
103 _RL dz1,dz2
104 _RL radius,StartingFlux
105 _RL dtts,dt
106 C Arrays
107 _RL Paa(Nr,nTracer)
108 _RL wda(Nr), mda(Nr), pda(Nr,nTracer)
109 C
110 C Pd, Wd :: tracers, vertical velocity in plume
111 C Md :: plume mass flux (?)
112 C Ad :: fractional area covered by plume
113 C Dd :: density in plume
114 C De :: density of environment
115 C PlumeEntrainment ::
116 _RL Ad(Nr),Wd(Nr),Dd(Nr),Md(Nr)
117 _RL De(Nr)
118 _RL PlumeEntrainment(Nr)
119 _RL Pd(Nr,nTracer)
120 CEOP
121
122 C-- Check to see if should convect now
123 C IF ( DIFFERENT_MULTIPLE(cAdjFreq,myTime,deltaTClock) ) THEN
124 IF ( .true. ) THEN
125 C local initialization
126
127 C Copy some arrays
128 dtts = dTtracerLev(1)
129
130 C start k-loop
131
132 DO k=1,KMax-1
133
134 C initialize the plume T,S,density, and w velocity
135
136 DO ktr=1,nTracerInUse
137 Pd(k,ktr) = tracerEnv(k,ktr)
138 ENDDO
139 Dd(k)=STATE1(Pd(k,2),Pd(k,1),i,j,k,bi,bj,myThid)
140 De(k)=Dd(k)
141 CML print *, 'ml-opps:', i,j,k,tracerEnv(k,2),tracerEnv(k,1),
142 CML & Dd(k),Pd(k,1),Pd(k,2)
143 CML compute vertical velocity at cell centers from GCM velocity
144 Wd(k)= - .5*(wVel(K)+wVel(K+1))
145 CML(
146 CML avoid division by zero
147 CML IF (Wd(K) .EQ. 0.D0) Wd(K) = 2.23e-16
148 CML)
149
150 C guess at initial top grid cell vertical velocity
151
152 CML Wd(k) = 0.03
153
154 C these estimates of initial plume velocity based on plume size and
155 C top grid cell water mass
156
157 c Wd(k) = 0.5*drF(k)/(dtts*FRACTIONAL_AREA)
158 c Wd(k) = 0.5*drF(k)/dtts
159
160 wsqr=Wd(k)*Wd(k)
161 PlumeEntrainment(k) = 0.0
162
163 #ifdef ALLOW_OPPS_DEBUG
164 IF ( OPPSdebugLevel.GE.debLevB ) THEN
165 WRITE(msgBuf,'(A,I3)')
166 & 'S/R OPPS_CALC: doing old lowerparcel', k
167 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
168 & SQUEEZE_RIGHT , 1)
169 ENDIF
170 #endif /* ALLOW_OPPS_DEBUG */
171 radius=PlumeRadius
172 StartingFlux=radius*radius*Wd(k)*Dd(k)
173 oldflux=StartingFlux
174
175 dz2=DrF(k)
176 DO k2=k,KMax-1
177 D1=STATE1( Pd(k2,2), Pd(k2,1),i,j,k2+1,bi,bj,myThid)
178 D2=STATE1( tracerEnv(k2+1,2), tracerEnv(k2+1,1),
179 & i,j,k2+1,bi,bj,myThid)
180 De(k2+1)=D2
181
182 C To start downward, parcel has to initially be heavier than environment
183 C but after it has started moving, we continue plume until plume tke or
184 C flux goes negative
185
186 CML & _hFacC(i,j,k-1,bi,bj)
187 CML & *_hFacC(i,j,k,bi,bj) .GT. 0.
188 CML & .AND.
189 IF (D2-D1 .LT. STABILITY_THRESHOLD.or.k2.ne.k) THEN
190 dz1=dz2
191 dz2=DrF(k2+1)
192
193 C find mass flux according to eq.(3) from paper by vertical integration
194
195 newflux=oldflux+e2*radius*Wd(k2)*Dd(k2)*
196 & .5*(dz1+dz2)
197 CML print *, 'ml-opps:', i,j,k,oldflux,newflux,e2,radius,
198 CML & Wd(k2),Dd(k2),Pd(k2,1),Pd(k2,2),dz1,dz2
199
200 PlumeEntrainment(k2+1) = newflux/StartingFlux
201
202 IF(newflux.LE.0.0) then
203 #ifdef ALLOW_OPPS_DEBUG
204 IF ( OPPSdebugLevel.GE.debLevA ) THEN
205 WRITE(msgBuf,'(A,I3)')
206 & 'S/R OPPS_CALC: Plume entrained to zero at level ', k2
207 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
208 & SQUEEZE_RIGHT , 1)
209 ENDIF
210 #endif /* ALLOW_OPPS_DEBUG */
211 maxdepth = k2
212 if(maxdepth.eq.k) goto 1000
213 goto 1
214 endif
215
216 C entrainment rate is basically a scaled mass flux dM/M
217
218 entrainrate = (newflux - oldflux)/newflux
219 oldflux = newflux
220
221 C mix var(s) are the average environmental values over the two grid levels
222
223 DO ktr=1,nTracerInUse
224 pmix=(dz1*tracerEnv(k2,ktr)+dz2*tracerEnv(k2+1,ktr))
225 & /(dz1+dz2)
226 Pd(k2+1,ktr)=Pd(k2,ktr)
227 & - entrainrate*(pmix - Pd(k2,ktr))
228 ENDDO
229
230 C compute the density at this level for the buoyancy term in the
231 C vertical k.e. equation
232
233 Dd(k2+1)=STATE1(Pd(k2+1,2),Pd(k2+1,1),i,j,k2+1,bi,bj,myThid)
234
235 C next, solve for the vertical velocity k.e. using combined eq. (4)
236 C and eq (5) from the paper
237
238 #ifdef ALLOW_OPPS_DEBUG
239 IF ( OPPSdebugLevel.GE.debLevA ) THEN
240 WRITE(msgBuf,'(A,3E12.4,I3)')
241 & 'S/R OPPS_CALC: Dd,De,entr,k ',Dd(k2),De(k2),entrainrate,k2
242 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
243 & SQUEEZE_RIGHT , 1)
244 ENDIF
245 #endif /* ALLOW_OPPS_DEBUG */
246 CML insert Eq. (4) into Eq. (5) to get something like this for wp^2
247 wsqr = wsqr - wsqr*abs(entrainrate)+ gravity*
248 & (dz1*(Dd(k2)-De(k2))/De(k2)
249 & +dz2*(Dd(k2+1)-De(k2+1))/De(k2+1))
250
251 C if negative k.e. then plume has reached max depth, get out of loop
252
253 IF(wsqr.LE.0.0)then
254 maxdepth = k2
255 #ifdef ALLOW_OPPS_DEBUG
256 IF ( OPPSdebugLevel.GE.debLevA ) THEN
257 WRITE(msgBuf,'(A,I3)')
258 & 'S/R OPPS_CALC: Plume velocity went to zero at level ', k2
259 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
260 & SQUEEZE_RIGHT , 1)
261 WRITE(msgBuf,'(A,4A14)')
262 & 'S/R OPPS_CALC: ', 'wsqr', 'entrainrate',
263 & '(Dd-De)/De up', '(Dd-De)/De do'
264 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
265 & SQUEEZE_RIGHT , 1)
266 WRITE(msgBuf,'(A,4E14.6)')
267 & 'S/R OPPS_CALC: ', wsqr, entrainrate,
268 & (Dd(k2)-De(k2))/De(k2), (Dd(k2+1)-De(k2+1))/De(k2+1)
269 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
270 & SQUEEZE_RIGHT , 1)
271 ENDIF
272 #endif /* ALLOW_OPPS_DEBUG */
273 if(maxdepth.eq.k) goto 1000
274 goto 1
275 endif
276 Wd(k2+1)=sqrt(wsqr)
277
278 C compute a new radius based on the new mass flux at this grid level
279 C from Eq. (4)
280
281 radius=sqrt(newflux/(Wd(k2)*Dd(k2)))
282 ELSE
283 maxdepth=k2
284 if(maxdepth.eq.k) goto 1000
285 GOTO 1
286 ENDIF
287 ENDDO
288
289 C plume has reached the bottom
290
291 MaxDepth=kMax
292
293 1 CONTINUE
294
295 Ad(k)=FRACTIONAL_AREA
296 IC=0
297
298 C start iteration on fractional area, not used in OGCM implementation
299
300 DO IC=1,Max_ABE_Iterations
301
302 C next compute the mass flux beteen each grid box using the entrainment
303
304 Md(k)=Wd(k)*Ad(k)
305
306 DO k2=k+1,maxDepth
307 Md(k2)=Md(k)*PlumeEntrainment(k2)
308 #ifdef ALLOW_OPPS_DEBUG
309 IF ( OPPSdebugLevel.GE.debLevA ) THEN
310 WRITE(msgBuf,'(A,2E12.4,I3)')
311 & 'S/R OPPS_CALC: Md, Wd, and k are ',Md(k2),Wd(k2),k2
312 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
313 & SQUEEZE_RIGHT , 1)
314 ENDIF
315 #endif /* ALLOW_OPPS_DEBUG */
316 ENDDO
317
318 C Now move on to calculate new temperature using flux from
319 C Td, Sd, Wd, ta, sa, and we. Values for these variables are at
320 C center of grid cell, use weighted average to get boundary values
321
322 C use a timestep limited by the GCM model timestep and the maximum plume
323 C velocity (CFL criteria)
324
325 C calculate the weighted wd, td, and sd
326 dt = dtts
327 do k2=k,maxDepth-1
328 IF ( Wd(K2) .NE. 0. _d 0 ) dt = min(dt,drF(k2)/Wd(k2))
329
330 C time integration will be integer number of steps to get one
331 C gcm time step
332
333 ntime = nint(0.5*int(dtts/dt))
334 if(ntime.eq.0) then
335 ntime = 1
336 endif
337
338 C make sure area weighted vertical velocities match; in other words
339 C make sure mass in equals mass out at the intersection of each grid
340 C cell. Eq. (20)
341
342 mda(k2) = (md(k2)*drF(k2)+md(k2+1)*drF(k2+1))/
343 * (drF(k2)+drF(k2+1))
344
345 wda(k2) = (wd(k2)*drF(k2)+wd(k2+1)*drF(k2+1))/
346 * (drF(k2)+drF(k2+1))
347
348 DO ktr = 1, nTracerInUse
349 Pda(k2,ktr) = Pd(k2,ktr)
350 Paa(k2,ktr) = tracerEnv(k2+1,ktr)
351 ENDDO
352
353 enddo
354 dt = min(dt,dtts)
355 #ifdef ALLOW_OPPS_DEBUG
356 IF ( OPPSdebugLevel.GE.debLevA ) THEN
357 WRITE(msgBuf,'(A,F14.4)')
358 & 'S/R OPPS_CALC: time step = ', dt
359 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
360 & SQUEEZE_RIGHT , 1)
361 ENDIF
362 #endif /* ALLOW_OPPS_DEBUG */
363 DO ktr=1,nTracerInUse
364 Pda(maxdepth,ktr) = Pd(maxdepth,ktr)
365 ENDDO
366
367 kmx = maxdepth-1
368 do nn=1,ntime
369
370 C top point
371
372 DO ktr = 1,nTracerInUse
373 tracerEnv(k,ktr) = tracerEnv(k,ktr)-
374 & (mda(k)*(Pda(k,ktr)-Paa(k,ktr)))*dt*recip_drF(k)
375 ENDDO
376
377 C now do inner points if there are any
378
379 CML if(Maxdepth-k.gt.1) then
380 CML This if statement is superfluous
381 CML IF ( k .LT. Maxdepth-1 ) THEN
382 CML DO k2=k+1,Maxdepth-1
383 CML mda(maxDepth) = 0.
384 DO k2=k+1,kmx
385 k2m1 = max(k,k2-1)
386 k2p1 = max(k2+1,maxDepth)
387
388 DO ktr = 1,nTracerInUse
389 tracerEnv(k2,ktr) = tracerEnv(k2,ktr) +
390 & (mda(k2m1)*(Pda(k2m1,ktr)-Paa(k2m1,ktr))
391 & -mda(k2) *(Pda(k2,ktr) -Paa(k2,ktr)) )
392 & *dt*recip_drF(k2)
393 ENDDO
394 ENDDO
395 CML This if statement is superfluous
396 CML ENDIF
397
398 C bottom point
399
400 DO ktr=1,nTracerInUse
401 tracerEnv(kmx+1,ktr) = tracerEnv(kmx+1,ktr)+
402 & mda(kmx)*(Pda(kmx,ktr)-Paa(kmx,ktr))*dt*recip_drF(kmx+1)
403 ENDDO
404
405 C set the environmental temp and salinity to equal new fields
406
407 DO ktr=1,nTracerInUse
408 DO k2=1,kmx
409 paa(k2,ktr) = tracerEnv(k2+1,ktr)
410 ENDDO
411 ENDDO
412
413 C end loop on number of time integration steps
414
415 enddo
416 ENDDO
417
418 C count convection event in this grid cell
419
420 OPPSconvectCount(I,J,K,bi,bj) =
421 & OPPSconvectCount(I,J,K,bi,bj) + 1. _d 0
422
423 C jump here if k = maxdepth or if level not unstable, go to next
424 C profile point
425
426 1000 continue
427
428 C end of k-loop
429
430 ENDDO
431
432 C-- End IF (DIFFERENT_MULTIPLE)
433 ENDIF
434
435 RETURN
436 END
437
438 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
439
440 _RL FUNCTION STATE1( sLoc, tLoc, i, j, kRef, bi, bj, myThid )
441
442 C !DESCRIPTION: \bv
443 C *===============================================================*
444 C | o SUBROUTINE STATE1
445 C | Calculates rho(S,T,p)
446 C | It is absolutely necessary to compute
447 C | the full rho and not sigma=rho-rhoConst, because
448 C | density is used as a scale factor for fluxes and velocities
449 C *===============================================================*
450 C \ev
451
452 C !USES:
453 IMPLICIT NONE
454 C == Global variables ==
455 #include "SIZE.h"
456 #include "EEPARAMS.h"
457 #include "PARAMS.h"
458 #include "GRID.h"
459 #include "DYNVARS.h"
460 #ifdef ALLOW_NONHYDROSTATIC
461 # include "NH_VARS.h"
462 #endif /* ALLOW_NONHYDROSTATIC */
463
464 C !INPUT/OUTPUT PARAMETERS:
465 C == Routine arguments ==
466 INTEGER i, j, kRef, bi, bj, myThid
467 _RL tLoc, sLoc
468
469 C !LOCAL VARIABLES:
470 C == Local variables ==
471 _RL rhoLoc
472 _RL pLoc
473
474 CMLC estimate pressure from depth at cell centers
475 CML mtoSI = gravity*rhoConst
476 CML pLoc = ABS(rC(kRef))*mtoSI
477
478 IF ( usingZCoords ) THEN
479 C in Z coordinates the pressure is rho0 * (hydrostatic) Potential
480 #ifdef ALLOW_NONHYDROSTATIC
481 IF ( selectP_inEOS_Zc.EQ.3 ) THEN
482 pLoc = rhoConst*( totPhiHyd(i,j,kRef,bi,bj)
483 & + phi_nh(i,j,kRef,bi,bj)
484 & + phiRef(2*kRef)
485 & )*maskC(i,j,kRef,bi,bj)
486 ELSEIF ( selectP_inEOS_Zc.EQ.2 ) THEN
487 #else /* ALLOW_NONHYDROSTATIC */
488 IF ( selectP_inEOS_Zc.EQ.2 ) THEN
489 #endif /* ALLOW_NONHYDROSTATIC */
490 C----------
491 C NOTE: For now, totPhiHyd only contains the Potential anomaly
492 C since PhiRef has not (yet) been added in S/R DIAGS_PHI_HYD
493 C----------
494 pLoc = rhoConst*( totPhiHyd(i,j,kRef,bi,bj)
495 & + phiRef(2*kRef)
496 & )*maskC(i,j,kRef,bi,bj)
497 c ELSEIF ( selectP_inEOS_Zc.EQ.1 ) THEN
498 C note: for the case selectP_inEOS_Zc=0, also use pRef4EOS (set to
499 C rhoConst*phiRef(2*k) ) to reproduce same previous machine truncation
500 ELSEIF ( selectP_inEOS_Zc.LE.1 ) THEN
501 pLoc = pRef4EOS(kRef)*maskC(i,j,kRef,bi,bj)
502 ELSE
503 pLoc = rhoConst*phiRef(2*kRef)*maskC(i,j,kRef,bi,bj)
504 ENDIF
505 ELSEIF ( usingPCoords ) THEN
506 C in P coordinates the pressure is just the coordinate of the tracer point
507 pLoc = rC(kRef)* maskC(i,j,kRef,bi,bj)
508 ENDIF
509
510 CALL FIND_RHO_SCALAR( tLoc, sLoc, pLoc, rhoLoc, myThid )
511 STATE1 = rhoLoc
512
513 #endif /* ALLOW_OPPS */
514 RETURN
515 END
516
517 #ifdef OPPS_ORGCODE
518 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
519
520 C Listed below is the subroutine for use in parallel 3-d circulation code.
521 C It has been used in the parallel semtner-chervin code and is now being used
522 C In the POP code. The subroutine is called nlopps (long story to explain why).
523
524 C I've attached the version of lopps that we've been using in the simulations.
525 C There is one common block that is different from the standard model commons
526 C (countc) and it is not needed if the array convadj is not used. The routine
527 C does need "kmp" which is why the boundc common is included. For massively
528 C parallel codes (like POP) we think this will work well when converted from a
529 C "slab" (i=is,ie) to a column, which just means removing the "do i=is,ie" loop.
530 C There are differences between this
531 C code and the 1-d code and the earlier scheme implemented in 3-d models. These
532 C differences are described below.
533
534 subroutine nlopps(j,is,ie,ta,sa,gcmdz)
535
536 parameter (imt = 361 , jmt = 301 , km = 30 )
537
538 C Nlopps: E. Skyllingstad and T. Paluszkiewicz
539
540 C Version: December 11, 1996
541
542 C Nlopps: This version of lopps is significantly different from
543 C the original code developed by R. Romea and T. Paluskiewicz. The
544 C code uses a flux constraint to control the change in T and S at
545 C each grid level. First, a plume profile of T,S, and W are
546 C determined using the standard plume model, but with a detraining
547 C mass instead of entraining. Thus, the T and S plume
548 C characteristics still change, but the plume contracts in size
549 C rather than expanding ala classical entraining plumes. This
550 C is heuristically more in line with large eddy simulation results.
551 C At each grid level, the convergence of plume velocity determines
552 C the flux of T and S, which is conserved by using an upstream
553 C advection. The vertical velocity is balanced so that the area
554 C weighted upward velocity equals the area weighted downdraft
555 C velocity, ensuring mass conservation. The present implementation
556 C adjusts the plume for a time period equal to the time for 1/2 of
557 C the mass of the fastest moving level to move downward. As a
558 C consequence, the model does not completely adjust the profile at
559 C each model time step, but provides a smooth adjustment over time.
560
561 c include "params.h"
562 c include "plume_fast_inc.h"
563 c include "plume_fast.h"
564 c #include "loppsd.h"
565
566 real ta(imt,km),sa(imt,km),gcmdz(km),dz(km)
567 real pdensity,wsqr,oldflux,newflux,entrainrate,adtemp
568 REAL Del,D,dza1,dza2,kd,kd1,Smix,Thmix,PlumeS,PlumeT,PlumeD
569
570 INTEGER i,j,k
571 Clfh
572 integer is,ie,k2
573 Clfh
574 REAL D1,D2,state1,Density
575 REAL dz1,dz2
576 REAL radius,StartingFlux
577 real ttemp(km),stemp(km),taa(km),saa(km)
578 real wda(km),tda(km),sda(km),mda(km)
579 real dtts,dt,sumo,sumn
580 integer ntime,nn,kmx,ic
581
582 LOGICAL debug,done
583 INTEGER MAX_ABE_ITERATIONS
584 PARAMETER(MAX_ABE_ITERATIONS=1)
585 REAL PlumeRadius
586 REAL STABILITY_THRESHOLD
587 REAL FRACTIONAL_AREA
588 REAL MAX_FRACTIONAL_AREA
589 REAL VERTICAL_VELOCITY
590 REAL ENTRAINMENT_RATE
591 REAL e2
592 PARAMETER ( PlumeRadius = 100.D0 )
593 PARAMETER ( STABILITY_THRESHOLD = -1.E-4 )
594 PARAMETER ( FRACTIONAL_AREA = .1E0 )
595 PARAMETER ( MAX_FRACTIONAL_AREA = .8E0 )
596 PARAMETER ( VERTICAL_VELOCITY = .02E0 )
597 PARAMETER ( ENTRAINMENT_RATE = -.05E0 )
598 PARAMETER ( e2 = 2.E0*ENTRAINMENT_RATE )
599 ! Arrays.
600 REAL Ad(km),Sd(km),Td(km),Wd(km),Dd(km),Md(km)
601 REAL Se(km),Te(km),We(km),De(km)
602 REAL PlumeEntrainment(km)
603 REAL GridThickness(km)
604
605 C input kmp through a common block
606
607 common / boundc / wsx(imt,jmt),wsy(imt,jmt),hfs(imt,jmt),
608 1 ple(imt,jmt),kmp(imt,jmt),kmq(imt,jmt)
609 cwmseas
610 & ,wsx1(imt,jmt),wsy1(imt,jmt)
611 1 ,wsx2(imt,jmt),wsy2(imt,jmt)
612
613 C input the variables through a common
614
615 logical problem
616 common /countc/ convadj(imt,jmt,km),ics,depth(km),problem
617
618 C-----may want to setup an option to get this only on first call
619 C otherwise it is repetive
620 C griddz is initialize by call to setupgrid
621
622 dtts = 2400
623
624 do k=1,km
625 dz(k) = 0.01*gcmdz(k)
626 enddo
627
628 do k=1,km
629 GridThickness(k) = dz(k)
630 enddo
631
632 C modified to loop over slab
633
634 DO i=is,ie
635
636 numgridpoints=kmp(i,j)
637
638 C go to next column if only 1 grid point or on land
639
640 if(numgridpoints.le.1) goto 1100
641
642 C loop over depth
643
644 debug = .false.
645
646 C first save copy of initial profile
647
648 DO k=1,NumGridPoints
649 stemp(k)=sa(i,k)
650 ttemp(k)=ta(i,k)
651
652 C do a check of t and s range, if out of bounds set flag
653
654 if(problem) then
655 write(0,*)"Code in trouble before this nlopps call"
656 return
657 endif
658
659 if(sa(i,k).gt.40..or.ta(i,k).lt.-4.0) then
660 problem = .true.
661 write(0,*)"t out of range at j ",j
662 debug = .true.
663 return
664 endif
665 ENDDO
666
667 if(debug) then
668 write(*,*)"T and S Profile at ",i,j
669 write(*,*)(ta(i,k),sa(i,k),k=1,NumGridPoints)
670 endif
671
672 DO k=1,NumGridPoints-1
673
674 C initialize the plume T,S,density, and w velocity
675
676 Sd(k)=stemp(k)
677 Td(k)=ttemp(k)
678 Dd(k)=state1(stemp(k),ttemp(k),k)
679 De(k)=Dd(k)
680 c Wd(k)=VERTICAL_VELOCITY
681
682 C guess at initial top grid cell vertical velocity
683
684 Wd(k) = 0.03
685
686 C these estimates of initial plume velocity based on plume size and
687 C top grid cell water mass
688
689 c Wd(k) = 0.5*dz(k)/(dtts*FRACTIONAL_AREA)
690 c Wd(k) = 0.5*dz(k)/dtts
691
692 wsqr=Wd(k)*Wd(k)
693 PlumeEntrainment(k) = 0.0
694
695 if(debug) write(0,*) 'Doing old lowerparcel'
696 radius=PlumeRadius
697 StartingFlux=radius*radius*Wd(k)*Dd(k)
698 oldflux=StartingFlux
699
700 dz2=GridThickness(k)
701 DO k2=k,NumGridPoints-1
702 D1=state1(Sd(k2),Td(k2),k2+1)
703 D2=state1(stemp(k2+1),ttemp(k2+1),k2+1)
704 De(k2+1)=D2
705
706 C To start downward, parcel has to initially be heavier than environment
707 C but after it has started moving, we continue plume until plume tke or
708 C flux goes negative
709
710 IF (D2-D1 .LT. STABILITY_THRESHOLD.or.k2.ne.k) THEN
711 dz1=dz2
712 dz2=GridThickness(k2+1)
713
714 C define mass flux according to eq. 4 from paper
715
716 newflux=oldflux+e2*radius*Wd(k2)*Dd(k2)*0.50*
717 . (dz1+dz2)
718
719 PlumeEntrainment(k2+1) = newflux/StartingFlux
720
721 IF(newflux.LT.0.0) then
722 if(debug) then
723 write(0,*)"Plume entrained to zero at ",k2
724 endif
725 maxdepth = k2
726 if(maxdepth.eq.k) goto 1000
727 goto 1
728 endif
729
730 C entrainment rate is basically a scaled mass flux dM/M
731
732 entrainrate = (newflux - oldflux)/newflux
733 oldflux = newflux
734
735 C mix var(s) are the average environmental values over the two grid levels
736
737 smix=(dz1*stemp(k2)+dz2*stemp(k2+1))/(dz1+dz2)
738 thmix=(dz1*ttemp(k2)+dz2*ttemp(k2+1))/(dz1+dz2)
739
740 C first compute the new salinity and temperature for this level
741 C using equations 3.6 and 3.7 from the paper
742
743 sd(k2+1)=sd(k2) - entrainrate*(smix - sd(k2))
744 td(k2+1)=td(k2) - entrainrate*(thmix - td(k2))
745
746 C compute the density at this level for the buoyancy term in the
747 C vertical k.e. equation
748
749 Dd(k2+1)=state1(Sd(k2+1),Td(k2+1),k2+1)
750
751 C next, solve for the vertical velocity k.e. using combined eq. 4
752 C and eq 5 from the paper
753
754 if(debug) then
755 write(0,*)"Dd,De,entr,k ",Dd(k2),De(k2),entrainrate,k2
756 endif
757
758 wsqr = wsqr - wsqr*abs(entrainrate)+ 9.81*
759 . (dz1*(Dd(k2)-De(k2))/De(k2)
760 . +dz2*(Dd(k2+1)-De(k2+1))/De(k2+1))
761
762 C if negative k.e. then plume has reached max depth, get out of loop
763
764 IF(wsqr.LT.0.0)then
765 maxdepth = k2
766 if(debug) then
767 write(0,*)"Plume velocity went to zero at ",k2
768 endif
769 if(maxdepth.eq.k) goto 1000
770 goto 1
771 endif
772 Wd(k2+1)=sqrt(wsqr)
773
774 C compute a new radius based on the new mass flux at this grid level
775
776 radius=sqrt(newflux/(Wd(k2)*Dd(k2)))
777 ELSE
778 maxdepth=k2
779 if(maxdepth.eq.k) goto 1000
780 GOTO 1
781 ENDIF
782 ENDDO
783
784 C plume has reached the bottom
785
786 MaxDepth=NumGridPoints
787
788 1 continue
789
790 Ad(k)=FRACTIONAL_AREA
791 IC=0
792
793 C start iteration on fractional area, not used in OGCM implementation
794
795 DO IC=1,Max_ABE_Iterations
796
797 C next compute the mass flux beteen each grid box using the entrainment
798
799 92 continue
800 Md(k)=Wd(k)*Ad(k)
801
802 DO k2=k+1,MaxDepth
803 Md(k2)=Md(k)*PlumeEntrainment(k2)
804 if(debug) then
805 write(0,*)"Md, Wd, and k are ",Md(k2),Wd(k2),k2
806 endif
807 ENDDO
808
809 C Now move on to calculate new temperature using flux from
810 C Td, Sd, Wd, ta, sa, and we. Values for these variables are at
811 C center of grid cell, use weighted average to get boundary values
812
813 C use a timestep limited by the GCM model timestep and the maximum plume
814 C velocity (CFL criteria)
815
816 C calculate the weighted wd, td, and sd
817
818 dt = dtts
819 do k2=k,maxdepth-1
820 dt = min(dt,dz(k2)/wd(k2))
821
822 C time integration will be integer number of steps to get one
823 C gcm time step
824
825 ntime = nint(0.5*int(dtts/dt))
826 if(ntime.eq.0) then
827 ntime = 1
828 endif
829
830 C make sure area weighted vertical velocities match; in other words
831 C make sure mass in equals mass out at the intersection of each grid
832 C cell.
833
834 mda(k2) = (md(k2)*dz(k2)+md(k2+1)*dz(k2+1))/
835 * (dz(k2)+dz(k2+1))
836
837 wda(k2) = (wd(k2)*dz(k2)+wd(k2+1)*dz(k2+1))/
838 * (dz(k2)+dz(k2+1))
839
840 tda(k2) = td(k2)
841 sda(k2) = sd(k2)
842
843 taa(k2) = ttemp(k2+1)
844 saa(k2) = stemp(k2+1)
845
846 enddo
847 dt = min(dt,dtts)
848 if(debug) then
849 write(0,*)"Time step is ", dt
850 endif
851 tda(maxdepth) = td(maxdepth)
852 sda(maxdepth) = sd(maxdepth)
853
854 C do top and bottom points first
855
856 kmx = maxdepth-1
857 do nn=1,ntime
858
859 ttemp(k) = ttemp(k)-
860 * (mda(k)*(tda(k)-taa(k)))*dt*recip_drF(k)
861
862 stemp(k) = stemp(k)-
863 * (mda(k)*(sda(k)-saa(k)))*dt*recip_drF(k)
864
865 C now do inner points if there are any
866
867 if(Maxdepth-k.gt.1) then
868 do k2=k+1,Maxdepth-1
869
870 ttemp(k2) = ttemp(k2) +
871 * (mda(k2-1)*(tda(k2-1)-taa(k2-1))-
872 * mda(k2)*(tda(k2)-taa(k2)))
873 * *dt*recip_drF(k2)
874
875 stemp(k2) = stemp(k2) +
876 * (mda(k2-1)*(sda(k2-1)-saa(k2-1))-
877 * mda(k2)*(sda(k2)-saa(k2)))
878 * *dt*recip_drF(k2)
879
880 enddo
881 endif
882 ttemp(kmx+1) = ttemp(kmx+1)+
883 * (mda(kmx)*(tda(kmx)-taa(kmx)))*
884 * dt*recip_drF(kmx+1)
885
886 stemp(kmx+1) = stemp(kmx+1)+
887 * (mda(kmx)*(sda(kmx)-saa(kmx)))*
888 * dt*recip_drF(kmx+1)
889
890 C set the environmental temp and salinity to equal new fields
891
892 do k2=1,maxdepth-1
893 taa(k2) = ttemp(k2+1)
894 saa(k2) = stemp(k2+1)
895 enddo
896
897 C end loop on number of time integration steps
898
899 enddo
900 ENDDO
901 999 continue
902
903 C assume that it converged, so update the ta and sa with new fields
904
905 c if(i.gt.180.and.j.gt.200.and.i.lt.240) then
906 c write(*,*)"Converged ",i,j,k,maxdepth,ttemp(k+1),ta(i,k+1)
907 c endif
908 do k2=k,maxdepth
909 convadj(i,j,k2) = convadj(i,j,k2) + (ttemp(k2)-
910 * ta(i,k2))
911 sa(i,k2) = stemp(k2)
912 ta(i,k2) = ttemp(k2)
913 c if(i.gt.180.and.j.gt.200.and.i.lt.240) then
914 c write(*,*)"convadj ",convadj(i,j,k2)
915 c endif
916
917 C see if nlopps messed up
918
919 if(sa(i,k).gt.40..or.ta(i,k).lt.-4.0) then
920 problem = .true.
921 write(0,*)"t out of range at j after adjust",j
922 debug = .true.
923 endif
924
925 enddo
926
927 C jump here if k = maxdepth or if level not unstable, go to next
928 C profile point
929
930 1000 continue
931
932 C end loop on k, move on to next possible plume
933
934 ENDDO
935 1100 continue
936
937 C i loop
938
939 ENDDO
940 END
941
942 #endif /* OPPS_ORGCODE */

  ViewVC Help
Powered by ViewVC 1.1.22