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