1 |
C $Header: /u/gcmpack/MITgcm/pkg/opps/opps_calc.F,v 1.10 2011/07/19 13:08:24 mlosch 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 orginially 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 |
999 continue |
418 |
|
419 |
C count convection event in this grid cell |
420 |
|
421 |
OPPSconvectCount(I,J,K,bi,bj) = |
422 |
& OPPSconvectCount(I,J,K,bi,bj) + 1. _d 0 |
423 |
|
424 |
C jump here if k = maxdepth or if level not unstable, go to next |
425 |
C profile point |
426 |
|
427 |
1000 continue |
428 |
|
429 |
C end of k-loop |
430 |
|
431 |
ENDDO |
432 |
|
433 |
C-- End IF (DIFFERENT_MULTIPLE) |
434 |
ENDIF |
435 |
|
436 |
RETURN |
437 |
END |
438 |
|
439 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
440 |
|
441 |
_RL FUNCTION STATE1( sLoc, tLoc, i, j, kRef, bi, bj, myThid ) |
442 |
|
443 |
C !DESCRIPTION: \bv |
444 |
C *===============================================================* |
445 |
C | o SUBROUTINE STATE1 |
446 |
C | Calculates rho(S,T,p) |
447 |
C | It is absolutely necessary to compute |
448 |
C | the full rho and not sigma=rho-rhoConst, because |
449 |
C | density is used as a scale factor for fluxes and velocities |
450 |
C *===============================================================* |
451 |
C \ev |
452 |
|
453 |
C !USES: |
454 |
IMPLICIT NONE |
455 |
C == Global variables == |
456 |
#include "SIZE.h" |
457 |
#include "EEPARAMS.h" |
458 |
#include "PARAMS.h" |
459 |
#include "GRID.h" |
460 |
#include "DYNVARS.h" |
461 |
|
462 |
C !INPUT/OUTPUT PARAMETERS: |
463 |
C == Routine arguments == |
464 |
INTEGER i, j, kRef, bi, bj, myThid |
465 |
_RL tLoc, sLoc |
466 |
|
467 |
C !LOCAL VARIABLES: |
468 |
C == Local variables == |
469 |
_RL rhoLoc |
470 |
_RL pLoc |
471 |
|
472 |
CMLC estimate pressure from depth at cell centers |
473 |
CML mtoSI = gravity*rhoConst |
474 |
CML pLoc = ABS(rC(kRef))*mtoSI |
475 |
|
476 |
IF ( usingZCoords ) THEN |
477 |
C in Z coordinates the pressure is rho0 * (hydrostatic) Potential |
478 |
IF ( useDynP_inEos_Zc ) THEN |
479 |
C---------- |
480 |
C NOTE: For now, totPhiHyd only contains the Potential anomaly |
481 |
C since PhiRef is not available for Atmos and has not (yet) |
482 |
C been added in S/R DIAGS_PHI_HYD |
483 |
C---------- |
484 |
pLoc = rhoConst*( totPhiHyd(i,j,kRef,bi,bj) |
485 |
& -rC(kRef)*gravity |
486 |
& )*maskC(i,j,kRef,bi,bj) |
487 |
ELSE |
488 |
pLoc = -rhoConst*rC(kRef)*gravity*maskC(i,j,kRef,bi,bj) |
489 |
ENDIF |
490 |
ELSEIF ( usingPCoords ) THEN |
491 |
C in P coordinates the pressure is just the coordinate of |
492 |
C the tracer point |
493 |
pLoc = rC(kRef)* maskC(i,j,kRef,bi,bj) |
494 |
ENDIF |
495 |
|
496 |
CALL FIND_RHO_SCALAR( tLoc, sLoc, pLoc, rhoLoc, myThid ) |
497 |
STATE1 = rhoLoc |
498 |
|
499 |
#endif /* ALLOW_OPPS */ |
500 |
RETURN |
501 |
END |
502 |
|
503 |
#ifdef OPPS_ORGCODE |
504 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
505 |
|
506 |
C Listed below is the subroutine for use in parallel 3-d circulation code. |
507 |
C It has been used in the parallel semtner-chervin code and is now being used |
508 |
C In the POP code. The subroutine is called nlopps (long story to explain why). |
509 |
|
510 |
C I've attached the version of lopps that we've been using in the simulations. |
511 |
C There is one common block that is different from the standard model commons |
512 |
C (countc) and it is not needed if the array convadj is not used. The routine |
513 |
C does need "kmp" which is why the boundc common is included. For massively |
514 |
C parallel codes (like POP) we think this will work well when converted from a |
515 |
C "slab" (i=is,ie) to a column, which just means removing the "do i=is,ie" loop. c There are differences between this |
516 |
C code and the 1-d code and the earlier scheme implemented in 3-d models. These c differences are described below. |
517 |
|
518 |
|
519 |
subroutine nlopps(j,is,ie,ta,sa,gcmdz) |
520 |
|
521 |
parameter (imt = 361 , jmt = 301 , km = 30 ) |
522 |
|
523 |
C Nlopps: E. Skyllingstad and T. Paluszkiewicz |
524 |
|
525 |
C Version: December 11, 1996 |
526 |
|
527 |
C Nlopps: This version of lopps is significantly different from |
528 |
C the original code developed by R. Romea and T. Paluskiewicz. The |
529 |
C code uses a flux constraint to control the change in T and S at |
530 |
C each grid level. First, a plume profile of T,S, and W are |
531 |
C determined using the standard plume model, but with a detraining |
532 |
C mass instead of entraining. Thus, the T and S plume |
533 |
C characteristics still change, but the plume contracts in size |
534 |
C rather than expanding ala classical entraining plumes. This |
535 |
C is heuristically more in line with large eddy simulation results. |
536 |
C At each grid level, the convergence of plume velocity determines |
537 |
C the flux of T and S, which is conserved by using an upstream |
538 |
C advection. The vertical velocity is balanced so that the area |
539 |
C weighted upward velocity equals the area weighted downdraft |
540 |
C velocity, ensuring mass conservation. The present implementation |
541 |
C adjusts the plume for a time period equal to the time for 1/2 of |
542 |
C the mass of the fastest moving level to move downward. As a |
543 |
C consequence, the model does not completely adjust the profile at |
544 |
C each model time step, but provides a smooth adjustment over time. |
545 |
|
546 |
c include "params.h" |
547 |
c include "plume_fast_inc.h" |
548 |
c include "plume_fast.h" |
549 |
c #include "loppsd.h" |
550 |
|
551 |
real ta(imt,km),sa(imt,km),gcmdz(km),dz(km) |
552 |
real pdensity,wsqr,oldflux,newflux,entrainrate,adtemp |
553 |
REAL Del,D,dza1,dza2,kd,kd1,Smix,Thmix,PlumeS,PlumeT,PlumeD |
554 |
|
555 |
INTEGER i,j,k |
556 |
Clfh |
557 |
integer is,ie,k2 |
558 |
Clfh |
559 |
REAL D1,D2,state1,Density |
560 |
REAL dz1,dz2 |
561 |
REAL radius,StartingFlux |
562 |
real ttemp(km),stemp(km),taa(km),saa(km) |
563 |
real wda(km),tda(km),sda(km),mda(km) |
564 |
real dtts,dt,sumo,sumn |
565 |
integer ntime,nn,kmx,ic |
566 |
|
567 |
LOGICAL debug,done |
568 |
INTEGER MAX_ABE_ITERATIONS |
569 |
PARAMETER(MAX_ABE_ITERATIONS=1) |
570 |
REAL PlumeRadius |
571 |
REAL STABILITY_THRESHOLD |
572 |
REAL FRACTIONAL_AREA |
573 |
REAL MAX_FRACTIONAL_AREA |
574 |
REAL VERTICAL_VELOCITY |
575 |
REAL ENTRAINMENT_RATE |
576 |
REAL e2 |
577 |
PARAMETER ( PlumeRadius = 100.D0 ) |
578 |
PARAMETER ( STABILITY_THRESHOLD = -1.E-4 ) |
579 |
PARAMETER ( FRACTIONAL_AREA = .1E0 ) |
580 |
PARAMETER ( MAX_FRACTIONAL_AREA = .8E0 ) |
581 |
PARAMETER ( VERTICAL_VELOCITY = .02E0 ) |
582 |
PARAMETER ( ENTRAINMENT_RATE = -.05E0 ) |
583 |
PARAMETER ( e2 = 2.E0*ENTRAINMENT_RATE ) |
584 |
! Arrays. |
585 |
REAL Ad(km),Sd(km),Td(km),Wd(km),Dd(km),Md(km) |
586 |
REAL Se(km),Te(km),We(km),De(km) |
587 |
REAL PlumeEntrainment(km) |
588 |
REAL GridThickness(km) |
589 |
|
590 |
C input kmp through a common block |
591 |
|
592 |
common / boundc / wsx(imt,jmt),wsy(imt,jmt),hfs(imt,jmt), |
593 |
1 ple(imt,jmt),kmp(imt,jmt),kmq(imt,jmt) |
594 |
cwmseas |
595 |
& ,wsx1(imt,jmt),wsy1(imt,jmt) |
596 |
1 ,wsx2(imt,jmt),wsy2(imt,jmt) |
597 |
|
598 |
C input the variables through a common |
599 |
|
600 |
logical problem |
601 |
common /countc/ convadj(imt,jmt,km),ics,depth(km),problem |
602 |
|
603 |
|
604 |
C-----may want to setup an option to get this only on first call |
605 |
C otherwise it is repetive |
606 |
C griddz is initialize by call to setupgrid |
607 |
|
608 |
dtts = 2400 |
609 |
|
610 |
do k=1,km |
611 |
dz(k) = 0.01*gcmdz(k) |
612 |
enddo |
613 |
|
614 |
do k=1,km |
615 |
GridThickness(k) = dz(k) |
616 |
enddo |
617 |
|
618 |
C modified to loop over slab |
619 |
|
620 |
DO i=is,ie |
621 |
|
622 |
numgridpoints=kmp(i,j) |
623 |
|
624 |
C go to next column if only 1 grid point or on land |
625 |
|
626 |
if(numgridpoints.le.1) goto 1100 |
627 |
|
628 |
C loop over depth |
629 |
|
630 |
debug = .false. |
631 |
|
632 |
C first save copy of initial profile |
633 |
|
634 |
DO k=1,NumGridPoints |
635 |
stemp(k)=sa(i,k) |
636 |
ttemp(k)=ta(i,k) |
637 |
|
638 |
C do a check of t and s range, if out of bounds set flag |
639 |
|
640 |
if(problem) then |
641 |
write(0,*)"Code in trouble before this nlopps call" |
642 |
return |
643 |
endif |
644 |
|
645 |
if(sa(i,k).gt.40..or.ta(i,k).lt.-4.0) then |
646 |
problem = .true. |
647 |
write(0,*)"t out of range at j ",j |
648 |
debug = .true. |
649 |
return |
650 |
endif |
651 |
ENDDO |
652 |
|
653 |
if(debug) then |
654 |
write(*,*)"T and S Profile at ",i,j |
655 |
write(*,*)(ta(i,k),sa(i,k),k=1,NumGridPoints) |
656 |
endif |
657 |
|
658 |
DO k=1,NumGridPoints-1 |
659 |
|
660 |
C initialize the plume T,S,density, and w velocity |
661 |
|
662 |
Sd(k)=stemp(k) |
663 |
Td(k)=ttemp(k) |
664 |
Dd(k)=state1(stemp(k),ttemp(k),k) |
665 |
De(k)=Dd(k) |
666 |
c Wd(k)=VERTICAL_VELOCITY |
667 |
|
668 |
C guess at initial top grid cell vertical velocity |
669 |
|
670 |
Wd(k) = 0.03 |
671 |
|
672 |
C these estimates of initial plume velocity based on plume size and |
673 |
C top grid cell water mass |
674 |
|
675 |
c Wd(k) = 0.5*dz(k)/(dtts*FRACTIONAL_AREA) |
676 |
c Wd(k) = 0.5*dz(k)/dtts |
677 |
|
678 |
wsqr=Wd(k)*Wd(k) |
679 |
PlumeEntrainment(k) = 0.0 |
680 |
|
681 |
if(debug) write(0,*) 'Doing old lowerparcel' |
682 |
radius=PlumeRadius |
683 |
StartingFlux=radius*radius*Wd(k)*Dd(k) |
684 |
oldflux=StartingFlux |
685 |
|
686 |
dz2=GridThickness(k) |
687 |
DO k2=k,NumGridPoints-1 |
688 |
D1=state1(Sd(k2),Td(k2),k2+1) |
689 |
D2=state1(stemp(k2+1),ttemp(k2+1),k2+1) |
690 |
De(k2+1)=D2 |
691 |
|
692 |
C To start downward, parcel has to initially be heavier than environment |
693 |
C but after it has started moving, we continue plume until plume tke or |
694 |
C flux goes negative |
695 |
|
696 |
IF (D2-D1 .LT. STABILITY_THRESHOLD.or.k2.ne.k) THEN |
697 |
dz1=dz2 |
698 |
dz2=GridThickness(k2+1) |
699 |
|
700 |
C define mass flux according to eq. 4 from paper |
701 |
|
702 |
newflux=oldflux+e2*radius*Wd(k2)*Dd(k2)*0.50* |
703 |
. (dz1+dz2) |
704 |
|
705 |
PlumeEntrainment(k2+1) = newflux/StartingFlux |
706 |
|
707 |
IF(newflux.LT.0.0) then |
708 |
if(debug) then |
709 |
write(0,*)"Plume entrained to zero at ",k2 |
710 |
endif |
711 |
maxdepth = k2 |
712 |
if(maxdepth.eq.k) goto 1000 |
713 |
goto 1 |
714 |
endif |
715 |
|
716 |
C entrainment rate is basically a scaled mass flux dM/M |
717 |
|
718 |
entrainrate = (newflux - oldflux)/newflux |
719 |
oldflux = newflux |
720 |
|
721 |
C mix var(s) are the average environmental values over the two grid levels |
722 |
|
723 |
smix=(dz1*stemp(k2)+dz2*stemp(k2+1))/(dz1+dz2) |
724 |
thmix=(dz1*ttemp(k2)+dz2*ttemp(k2+1))/(dz1+dz2) |
725 |
|
726 |
C first compute the new salinity and temperature for this level |
727 |
C using equations 3.6 and 3.7 from the paper |
728 |
|
729 |
sd(k2+1)=sd(k2) - entrainrate*(smix - sd(k2)) |
730 |
td(k2+1)=td(k2) - entrainrate*(thmix - td(k2)) |
731 |
|
732 |
C compute the density at this level for the buoyancy term in the |
733 |
C vertical k.e. equation |
734 |
|
735 |
Dd(k2+1)=state1(Sd(k2+1),Td(k2+1),k2+1) |
736 |
|
737 |
C next, solve for the vertical velocity k.e. using combined eq. 4 |
738 |
C and eq 5 from the paper |
739 |
|
740 |
if(debug) then |
741 |
write(0,*)"Dd,De,entr,k ",Dd(k2),De(k2),entrainrate,k2 |
742 |
endif |
743 |
|
744 |
wsqr = wsqr - wsqr*abs(entrainrate)+ 9.81* |
745 |
. (dz1*(Dd(k2)-De(k2))/De(k2) |
746 |
. +dz2*(Dd(k2+1)-De(k2+1))/De(k2+1)) |
747 |
|
748 |
C if negative k.e. then plume has reached max depth, get out of loop |
749 |
|
750 |
IF(wsqr.LT.0.0)then |
751 |
maxdepth = k2 |
752 |
if(debug) then |
753 |
write(0,*)"Plume velocity went to zero at ",k2 |
754 |
endif |
755 |
if(maxdepth.eq.k) goto 1000 |
756 |
goto 1 |
757 |
endif |
758 |
Wd(k2+1)=sqrt(wsqr) |
759 |
|
760 |
C compute a new radius based on the new mass flux at this grid level |
761 |
|
762 |
radius=sqrt(newflux/(Wd(k2)*Dd(k2))) |
763 |
ELSE |
764 |
maxdepth=k2 |
765 |
if(maxdepth.eq.k) goto 1000 |
766 |
GOTO 1 |
767 |
ENDIF |
768 |
ENDDO |
769 |
|
770 |
C plume has reached the bottom |
771 |
|
772 |
MaxDepth=NumGridPoints |
773 |
|
774 |
1 continue |
775 |
|
776 |
Ad(k)=FRACTIONAL_AREA |
777 |
IC=0 |
778 |
|
779 |
C start iteration on fractional area, not used in OGCM implementation |
780 |
|
781 |
DO IC=1,Max_ABE_Iterations |
782 |
|
783 |
C next compute the mass flux beteen each grid box using the entrainment |
784 |
|
785 |
92 continue |
786 |
Md(k)=Wd(k)*Ad(k) |
787 |
|
788 |
DO k2=k+1,MaxDepth |
789 |
Md(k2)=Md(k)*PlumeEntrainment(k2) |
790 |
if(debug) then |
791 |
write(0,*)"Md, Wd, and k are ",Md(k2),Wd(k2),k2 |
792 |
endif |
793 |
ENDDO |
794 |
|
795 |
C Now move on to calculate new temperature using flux from |
796 |
C Td, Sd, Wd, ta, sa, and we. Values for these variables are at |
797 |
C center of grid cell, use weighted average to get boundary values |
798 |
|
799 |
C use a timestep limited by the GCM model timestep and the maximum plume |
800 |
C velocity (CFL criteria) |
801 |
|
802 |
C calculate the weighted wd, td, and sd |
803 |
|
804 |
dt = dtts |
805 |
do k2=k,maxdepth-1 |
806 |
dt = min(dt,dz(k2)/wd(k2)) |
807 |
|
808 |
C time integration will be integer number of steps to get one |
809 |
C gcm time step |
810 |
|
811 |
ntime = nint(0.5*int(dtts/dt)) |
812 |
if(ntime.eq.0) then |
813 |
ntime = 1 |
814 |
endif |
815 |
|
816 |
C make sure area weighted vertical velocities match; in other words |
817 |
C make sure mass in equals mass out at the intersection of each grid |
818 |
C cell. |
819 |
|
820 |
mda(k2) = (md(k2)*dz(k2)+md(k2+1)*dz(k2+1))/ |
821 |
* (dz(k2)+dz(k2+1)) |
822 |
|
823 |
wda(k2) = (wd(k2)*dz(k2)+wd(k2+1)*dz(k2+1))/ |
824 |
* (dz(k2)+dz(k2+1)) |
825 |
|
826 |
tda(k2) = td(k2) |
827 |
sda(k2) = sd(k2) |
828 |
|
829 |
taa(k2) = ttemp(k2+1) |
830 |
saa(k2) = stemp(k2+1) |
831 |
|
832 |
enddo |
833 |
dt = min(dt,dtts) |
834 |
if(debug) then |
835 |
write(0,*)"Time step is ", dt |
836 |
endif |
837 |
tda(maxdepth) = td(maxdepth) |
838 |
sda(maxdepth) = sd(maxdepth) |
839 |
|
840 |
C do top and bottom points first |
841 |
|
842 |
kmx = maxdepth-1 |
843 |
do nn=1,ntime |
844 |
|
845 |
ttemp(k) = ttemp(k)- |
846 |
* (mda(k)*(tda(k)-taa(k)))*dt*recip_drF(k) |
847 |
|
848 |
stemp(k) = stemp(k)- |
849 |
* (mda(k)*(sda(k)-saa(k)))*dt*recip_drF(k) |
850 |
|
851 |
C now do inner points if there are any |
852 |
|
853 |
if(Maxdepth-k.gt.1) then |
854 |
do k2=k+1,Maxdepth-1 |
855 |
|
856 |
ttemp(k2) = ttemp(k2) + |
857 |
* (mda(k2-1)*(tda(k2-1)-taa(k2-1))- |
858 |
* mda(k2)*(tda(k2)-taa(k2))) |
859 |
* *dt*recip_drF(k2) |
860 |
|
861 |
stemp(k2) = stemp(k2) + |
862 |
* (mda(k2-1)*(sda(k2-1)-saa(k2-1))- |
863 |
* mda(k2)*(sda(k2)-saa(k2))) |
864 |
* *dt*recip_drF(k2) |
865 |
|
866 |
enddo |
867 |
endif |
868 |
ttemp(kmx+1) = ttemp(kmx+1)+ |
869 |
* (mda(kmx)*(tda(kmx)-taa(kmx)))* |
870 |
* dt*recip_drF(kmx+1) |
871 |
|
872 |
stemp(kmx+1) = stemp(kmx+1)+ |
873 |
* (mda(kmx)*(sda(kmx)-saa(kmx)))* |
874 |
* dt*recip_drF(kmx+1) |
875 |
|
876 |
C set the environmental temp and salinity to equal new fields |
877 |
|
878 |
do k2=1,maxdepth-1 |
879 |
taa(k2) = ttemp(k2+1) |
880 |
saa(k2) = stemp(k2+1) |
881 |
enddo |
882 |
|
883 |
C end loop on number of time integration steps |
884 |
|
885 |
enddo |
886 |
ENDDO |
887 |
999 continue |
888 |
|
889 |
C assume that it converged, so update the ta and sa with new fields |
890 |
|
891 |
c if(i.gt.180.and.j.gt.200.and.i.lt.240) then |
892 |
c write(*,*)"Converged ",i,j,k,maxdepth,ttemp(k+1),ta(i,k+1) |
893 |
c endif |
894 |
do k2=k,maxdepth |
895 |
convadj(i,j,k2) = convadj(i,j,k2) + (ttemp(k2)- |
896 |
* ta(i,k2)) |
897 |
sa(i,k2) = stemp(k2) |
898 |
ta(i,k2) = ttemp(k2) |
899 |
c if(i.gt.180.and.j.gt.200.and.i.lt.240) then |
900 |
c write(*,*)"convadj ",convadj(i,j,k2) |
901 |
c endif |
902 |
|
903 |
C see if nlopps messed up |
904 |
|
905 |
if(sa(i,k).gt.40..or.ta(i,k).lt.-4.0) then |
906 |
problem = .true. |
907 |
write(0,*)"t out of range at j after adjust",j |
908 |
debug = .true. |
909 |
endif |
910 |
|
911 |
enddo |
912 |
|
913 |
C jump here if k = maxdepth or if level not unstable, go to next |
914 |
C profile point |
915 |
|
916 |
1000 continue |
917 |
|
918 |
C end loop on k, move on to next possible plume |
919 |
|
920 |
ENDDO |
921 |
1100 continue |
922 |
|
923 |
C i loop |
924 |
|
925 |
ENDDO |
926 |
END |
927 |
|
928 |
#endif /* OPPS_ORGCODE */ |