1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "BLING_OPTIONS.h" |
5 |
|
6 |
C-- File bling_carbon_chem.F: |
7 |
C-- Contents |
8 |
C-- o CALC_PCO2 |
9 |
C-- o CALC_PCO2_APPROX |
10 |
C-- o CARBON_COEFFS |
11 |
C-- o CARBON_COEFFS_PRESSURE_DEP |
12 |
|
13 |
CBOP |
14 |
subroutine CALC_PCO2( |
15 |
I donewt,inewtonmax,ibrackmax, |
16 |
I t,s,diclocal,pt,sit,ta, |
17 |
I k1local,k2local, |
18 |
I k1plocal,k2plocal,k3plocal, |
19 |
I kslocal,kblocal,kwlocal, |
20 |
I ksilocal,kflocal, |
21 |
I k0local, fugflocal, |
22 |
I fflocal,btlocal,stlocal,ftlocal, |
23 |
U pHlocal,pCO2surfloc, |
24 |
I i,j,k,bi,bj,myIter,myThid ) |
25 |
|
26 |
C !DESCRIPTION: |
27 |
C surface ocean inorganic carbon chemistry to OCMIP2 |
28 |
C regulations modified from OCMIP2 code; |
29 |
C Mick Follows, MIT, Oct 1999. |
30 |
C Apr 2011: fix vapour bug (following Bennington) |
31 |
|
32 |
implicit none |
33 |
|
34 |
C === Global variables === |
35 |
#include "SIZE.h" |
36 |
#include "DYNVARS.h" |
37 |
#include "EEPARAMS.h" |
38 |
#include "PARAMS.h" |
39 |
#include "GRID.h" |
40 |
#include "FFIELDS.h" |
41 |
#include "BLING_VARS.h" |
42 |
|
43 |
C == Routine arguments == |
44 |
C diclocal = total inorganic carbon (mol/m^3) |
45 |
C where 1 T = 1 metric ton = 1000 kg |
46 |
C ta = total alkalinity (eq/m^3) |
47 |
C pt = inorganic phosphate (mol/^3) |
48 |
C sit = inorganic silicate (mol/^3) |
49 |
C t = temperature (degrees C) |
50 |
C s = salinity (PSU) |
51 |
INTEGER donewt |
52 |
INTEGER inewtonmax |
53 |
INTEGER ibrackmax |
54 |
_RL t, s, pt, sit, ta |
55 |
_RL pCO2surfloc, diclocal, pHlocal |
56 |
_RL fflocal, btlocal, stlocal, ftlocal |
57 |
_RL k1local, k2local |
58 |
_RL k1plocal, k2plocal, k3plocal |
59 |
_RL kslocal, kblocal, kwlocal, ksilocal, kflocal |
60 |
_RL k0local, fugflocal |
61 |
INTEGER i,j,k,bi,bj,myIter |
62 |
INTEGER myThid |
63 |
CEOP |
64 |
|
65 |
C == Local variables == |
66 |
C INPUT |
67 |
C phlo= lower limit of pH range |
68 |
C phhi= upper limit of pH range |
69 |
C atmpres = atmospheric pressure in atmospheres (1 atm==1013.25mbar) |
70 |
C OUTPUT |
71 |
C co2star = CO2*water (mol/m^3) |
72 |
C pco2surf = oceanic pCO2 (ppmv) |
73 |
C --------------------------------------------------------------------- |
74 |
C OCMIP NOTE: Some words about units - (JCO, 4/4/1999) |
75 |
C - Models carry tracers in mol/m^3 (on a per volume basis) |
76 |
C - Conversely, this routine, which was written by |
77 |
C observationalists (C. Sabine and R. Key), passes input |
78 |
C arguments in umol/kg (i.e., on a per mass basis) |
79 |
C - I have changed things slightly so that input arguments are in |
80 |
C mol/m^3, |
81 |
C - Thus, all input concentrations (diclocal, ta, pt, and st) should be |
82 |
C given in mol/m^3; output arguments "co2star" and "dco2star" |
83 |
C are likewise be in mol/m^3. |
84 |
C --------------------------------------------------------------------- |
85 |
|
86 |
_RL phhi |
87 |
_RL phlo |
88 |
_RL c |
89 |
_RL a |
90 |
_RL a2 |
91 |
_RL da |
92 |
_RL b |
93 |
_RL b2 |
94 |
_RL db |
95 |
_RL fn |
96 |
_RL df |
97 |
_RL deltax |
98 |
_RL x |
99 |
_RL x2 |
100 |
_RL x3 |
101 |
_RL xmid |
102 |
_RL ftest |
103 |
_RL htotal |
104 |
_RL htotal2 |
105 |
_RL co2star |
106 |
_RL phguess |
107 |
_RL fco2 |
108 |
INTEGER inewton |
109 |
INTEGER ibrack |
110 |
INTEGER hstep |
111 |
_RL fni(3) |
112 |
_RL xlo |
113 |
_RL xhi |
114 |
_RL xguess |
115 |
_RL k123p |
116 |
_RL k12p |
117 |
_RL k12 |
118 |
c --------------------------------------------------------------------- |
119 |
c import donewt flag |
120 |
c set donewt = 1 for newton-raphson iteration |
121 |
c set donewt = 0 for bracket and bisection |
122 |
c --------------------------------------------------------------------- |
123 |
C Change units from the input of mol/m^3 -> mol/kg: |
124 |
c (1 mol/m^3) x (1 m^3/1024.5 kg) |
125 |
c where the ocean mean surface density is 1024.5 kg/m^3 |
126 |
c Note: mol/kg are actually what the body of this routine uses |
127 |
c for calculations. Units are reconverted back to mol/m^3 at the |
128 |
c end of this routine. |
129 |
c --------------------------------------------------------------------- |
130 |
c To convert input in mol/m^3 -> mol/kg |
131 |
pt=pt*permil |
132 |
sit=sit*permil |
133 |
ta=ta*permil |
134 |
diclocal=diclocal*permil |
135 |
c --------------------------------------------------------------------- |
136 |
c set first guess and brackets for [H+] solvers |
137 |
c first guess (for newton-raphson) |
138 |
phguess = phlocal |
139 |
|
140 |
|
141 |
c bracketing values (for bracket/bisection) |
142 |
phhi = 10.0 |
143 |
phlo = 5.0 |
144 |
c convert to [H+]... |
145 |
xguess = 10.0**(-phguess) |
146 |
xlo = 10.0**(-phhi) |
147 |
xhi = 10.0**(-phlo) |
148 |
xmid = (xlo + xhi)*0.5 |
149 |
|
150 |
|
151 |
c---------------------------------------------------------------- |
152 |
c iteratively solve for [H+] |
153 |
c (i) Newton-Raphson method with fixed number of iterations, |
154 |
c use previous [H+] as first guess |
155 |
|
156 |
c select newton-raphson, inewt=1 |
157 |
c else select bracket and bisection |
158 |
|
159 |
cQQQQQ |
160 |
if( donewt .eq. 1)then |
161 |
c......................................................... |
162 |
c NEWTON-RAPHSON METHOD |
163 |
c......................................................... |
164 |
x = xguess |
165 |
cdiags |
166 |
c WRITE(0,*)'xguess ',xguess |
167 |
cdiags |
168 |
do inewton = 1, inewtonmax |
169 |
c set some common combinations of parameters used in |
170 |
c the iterative [H+] solvers |
171 |
x2=x*x |
172 |
x3=x2*x |
173 |
k12 = k1local*k2local |
174 |
k12p = k1plocal*k2plocal |
175 |
k123p = k12p*k3plocal |
176 |
c = 1.0 + stlocal/kslocal |
177 |
a = x3 + k1plocal*x2 + k12p*x + k123p |
178 |
a2=a*a |
179 |
da = 3.0*x2 + 2.0*k1plocal*x + k12p |
180 |
b = x2 + k1local*x + k12 |
181 |
b2=b*b |
182 |
db = 2.0*x + k1local |
183 |
|
184 |
c Evaluate f([H+]) and f_prime([H+]) |
185 |
c fn = hco3+co3+borate+oh+hpo4+2*po4+silicate+hfree |
186 |
c +hso4+hf+h3po4-ta |
187 |
fn = k1local*x*diclocal/b + |
188 |
& 2.0*diclocal*k12/b + |
189 |
& btlocal/(1.0 + x/kblocal) + |
190 |
& kwlocal/x + |
191 |
& pt*k12p*x/a + |
192 |
& 2.0*pt*k123p/a + |
193 |
& sit/(1.0 + x/ksilocal) - |
194 |
& x/c - |
195 |
& stlocal/(1.0 + kslocal/x/c) - |
196 |
& ftlocal/(1.0 + kflocal/x) - |
197 |
& pt*x3/a - |
198 |
& ta |
199 |
|
200 |
c df = dfn/dx |
201 |
cdiags |
202 |
c WRITE(0,*)'values',b2,kblocal,x2,a2,c,x |
203 |
cdiags |
204 |
df = ((k1local*diclocal*b) - k1local*x*diclocal*db)/b2 - |
205 |
& 2.0*diclocal*k12*db/b2 - |
206 |
& btlocal/kblocal/(1.0+x/kblocal)**2. - |
207 |
& kwlocal/x2 + |
208 |
& (pt*k12p*(a - x*da))/a2 - |
209 |
& 2.0*pt*k123p*da/a2 - |
210 |
& sit/ksilocal/(1.0+x/ksilocal)**2. + |
211 |
& 1.0/c + |
212 |
& stlocal*(1.0 + kslocal/x/c)**(-2.0)*(kslocal/c/x2) + |
213 |
& ftlocal*(1.0 + kflocal/x)**(-2.)*kflocal/x2 - |
214 |
& pt*x2*(3.0*a-x*da)/a2 |
215 |
c evaluate increment in [H+] |
216 |
deltax = - fn/df |
217 |
c update estimate of [H+] |
218 |
x = x + deltax |
219 |
cdiags |
220 |
c write value of x to check convergence.... |
221 |
c write(0,*)'inewton, x, deltax ',inewton, x, deltax |
222 |
c write(6,*) |
223 |
cdiags |
224 |
|
225 |
end do |
226 |
c end of newton-raphson method |
227 |
c.................................................... |
228 |
else |
229 |
c.................................................... |
230 |
C BRACKET AND BISECTION METHOD |
231 |
c.................................................... |
232 |
c (ii) If first step use Bracket and Bisection method |
233 |
c with fixed, large number of iterations |
234 |
do ibrack = 1, ibrackmax |
235 |
do hstep = 1,3 |
236 |
if(hstep .eq. 1)x = xhi |
237 |
if(hstep .eq. 2)x = xlo |
238 |
if(hstep .eq. 3)x = xmid |
239 |
c set some common combinations of parameters used in |
240 |
c the iterative [H+] solvers |
241 |
|
242 |
|
243 |
x2=x*x |
244 |
x3=x2*x |
245 |
k12 = k1local*k2local |
246 |
k12p = k1plocal*k2plocal |
247 |
k123p = k12p*k3plocal |
248 |
c = 1.0 + stlocal/kslocal |
249 |
a = x3 + k1plocal*x2 + k12p*x + k123p |
250 |
a2=a*a |
251 |
da = 3.0*x2 + 2.0*k1plocal*x + k12p |
252 |
b = x2 + k1local*x + k12 |
253 |
b2=b*b |
254 |
db = 2.0*x + k1local |
255 |
c evaluate f([H+]) for bracketing and mid-value cases |
256 |
fn = k1local*x*diclocal/b + |
257 |
& 2.0*diclocal*k12/b + |
258 |
& btlocal/(1.0 + x/kblocal) + |
259 |
& kwlocal/x + |
260 |
& pt*k12p*x/a + |
261 |
& 2.0*pt*k123p/a + |
262 |
& sit/(1.0 + x/ksilocal) - |
263 |
& x/c - |
264 |
& stlocal/(1.0 + kslocal/x/c) - |
265 |
& ftlocal/(1.0 + kflocal/x) - |
266 |
& pt*x3/a - |
267 |
& ta |
268 |
fni(hstep) = fn |
269 |
end do |
270 |
c now bracket solution within two of three |
271 |
ftest = fni(1)/fni(3) |
272 |
if(ftest .gt. 0.0)then |
273 |
xhi = xmid |
274 |
else |
275 |
xlo = xmid |
276 |
end if |
277 |
xmid = (xlo + xhi)*0.5 |
278 |
|
279 |
cdiags |
280 |
c write value of x to check convergence.... |
281 |
c WRITE(0,*)'bracket-bisection iteration ',ibrack, xmid |
282 |
cdiags |
283 |
end do |
284 |
c last iteration gives value |
285 |
x = xmid |
286 |
c end of bracket and bisection method |
287 |
c.................................... |
288 |
end if |
289 |
c iterative [H+] solver finished |
290 |
c---------------------------------------------------------------- |
291 |
|
292 |
c now determine pCO2 etc... |
293 |
c htotal = [H+], hydrogen ion conc |
294 |
htotal = x |
295 |
C Calculate [CO2*] as defined in DOE Methods Handbook 1994 Ver.2, |
296 |
C ORNL/CDIAC-74, dickson and Goyet, eds. (Ch 2 p 10, Eq A.49) |
297 |
htotal2=htotal*htotal |
298 |
co2star=diclocal*htotal2/(htotal2 + k1local*htotal |
299 |
& + k1local*k2local) |
300 |
phlocal=-log10(htotal) |
301 |
|
302 |
c --------------------------------------------------------------- |
303 |
c Add two output arguments for storing pCO2surf |
304 |
c Should we be using K0 or ff for the solubility here? |
305 |
c --------------------------------------------------------------- |
306 |
fco2 = co2star / k0local |
307 |
pCO2surfloc = fco2/fugflocal |
308 |
|
309 |
C ---------------------------------------------------------------- |
310 |
C Reconvert units back to original values for input arguments |
311 |
C no longer necessary???? |
312 |
C ---------------------------------------------------------------- |
313 |
c Reconvert from mol/kg -> mol/m^3 |
314 |
pt=pt/permil |
315 |
sit=sit/permil |
316 |
ta=ta/permil |
317 |
diclocal=diclocal/permil |
318 |
|
319 |
RETURN |
320 |
END |
321 |
|
322 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
323 |
|
324 |
CBOP |
325 |
subroutine CALC_PCO2_APPROX( |
326 |
I t,s,diclocal,pt,sit,ta, |
327 |
I k1local,k2local, |
328 |
I k1plocal,k2plocal,k3plocal, |
329 |
I kslocal,kblocal,kwlocal, |
330 |
I ksilocal,kflocal, |
331 |
I k0local, fugflocal, |
332 |
I fflocal,btlocal,stlocal,ftlocal, |
333 |
U pHlocal,pCO2surfloc,co3local, |
334 |
I i,j,k,bi,bj,myIter,myThid ) |
335 |
|
336 |
C !DESCRIPTION: |
337 |
C *==========================================================* |
338 |
C | SUBROUTINE CALC_PCO2_APPROX | |
339 |
C *==========================================================* |
340 |
C CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
341 |
C C New efficient pCO2 solver, Mick Follows CC |
342 |
C C Taka Ito CC |
343 |
C C Stephanie Dutkiewicz CC |
344 |
C C 20 April 2003 CC |
345 |
C CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
346 |
C Apr 2011: fix vapour bug (following Bennington) |
347 |
C Oct 2011: add CO3 extimation and pass out |
348 |
|
349 |
implicit none |
350 |
|
351 |
C === Global variables === |
352 |
#include "SIZE.h" |
353 |
#include "DYNVARS.h" |
354 |
#include "EEPARAMS.h" |
355 |
#include "PARAMS.h" |
356 |
#include "GRID.h" |
357 |
#include "FFIELDS.h" |
358 |
#include "BLING_VARS.h" |
359 |
|
360 |
C == Routine arguments == |
361 |
C diclocal = total inorganic carbon (mol/m^3) |
362 |
C where 1 T = 1 metric ton = 1000 kg |
363 |
C ta = total alkalinity (eq/m^3) |
364 |
C pt = inorganic phosphate (mol/^3) |
365 |
C sit = inorganic silicate (mol/^3) |
366 |
C t = temperature (degrees C) |
367 |
C s = salinity (PSU) |
368 |
_RL t, s, pt, sit, ta |
369 |
_RL pCO2surfloc, diclocal, pHlocal |
370 |
_RL fflocal, btlocal, stlocal, ftlocal |
371 |
_RL k1local, k2local |
372 |
_RL k1plocal, k2plocal, k3plocal |
373 |
_RL kslocal, kblocal, kwlocal, ksilocal, kflocal |
374 |
_RL k0local, fugflocal |
375 |
_RL co3local |
376 |
INTEGER i,j,k,bi,bj,myIter |
377 |
INTEGER myThid |
378 |
CEOP |
379 |
|
380 |
C == Local variables == |
381 |
_RL phguess |
382 |
_RL cag |
383 |
_RL bohg |
384 |
_RL hguess |
385 |
_RL stuff |
386 |
_RL gamm |
387 |
_RL hnew |
388 |
_RL co2s |
389 |
_RL h3po4g, h2po4g, hpo4g, po4g |
390 |
_RL siooh3g |
391 |
_RL fco2 |
392 |
|
393 |
|
394 |
c --------------------------------------------------------------------- |
395 |
C Change units from the input of mol/m^3 -> mol/kg: |
396 |
c (1 mol/m^3) x (1 m^3/1024.5 kg) |
397 |
c where the ocean mean surface density is 1024.5 kg/m^3 |
398 |
c Note: mol/kg are actually what the body of this routine uses |
399 |
c for calculations. Units are reconverted back to mol/m^3 at the |
400 |
c end of this routine. |
401 |
c To convert input in mol/m^3 -> mol/kg |
402 |
pt=pt*permil |
403 |
sit=sit*permil |
404 |
ta=ta*permil |
405 |
diclocal=diclocal*permil |
406 |
c --------------------------------------------------------------------- |
407 |
c set first guess and brackets for [H+] solvers |
408 |
c first guess (for newton-raphson) |
409 |
phguess = phlocal |
410 |
cmick - new approx method |
411 |
cmick - make estimate of htotal (hydrogen ion conc) using |
412 |
cmick appromate estimate of CA, carbonate alkalinity |
413 |
hguess = 10.0**(-phguess) |
414 |
cmick - first estimate borate contribution using guess for [H+] |
415 |
bohg = btlocal*kblocal/(hguess+kblocal) |
416 |
|
417 |
cmick - first estimate of contribution from phosphate |
418 |
cmick based on Dickson and Goyet |
419 |
stuff = hguess*hguess*hguess |
420 |
& + (k1plocal*hguess*hguess) |
421 |
& + (k1plocal*k2plocal*hguess) |
422 |
& + (k1plocal*k2plocal*k3plocal) |
423 |
h3po4g = (pt*hguess*hguess*hguess) / stuff |
424 |
h2po4g = (pt*k1plocal*hguess*hguess) / stuff |
425 |
hpo4g = (pt*k1plocal*k2plocal*hguess) / stuff |
426 |
po4g = (pt*k1plocal*k2plocal*k3plocal) / stuff |
427 |
|
428 |
cmick - estimate contribution from silicate |
429 |
cmick based on Dickson and Goyet |
430 |
siooh3g = sit*ksilocal / (ksilocal + hguess) |
431 |
|
432 |
cmick - now estimate carbonate alkalinity |
433 |
cag = ta - bohg - (kwlocal/hguess) + hguess |
434 |
& - hpo4g - 2.0 _d 0*po4g + h3po4g |
435 |
& - siooh3g |
436 |
|
437 |
cmick - now evaluate better guess of hydrogen ion conc |
438 |
cmick htotal = [H+], hydrogen ion conc |
439 |
gamm = diclocal/cag |
440 |
stuff = (1.0 _d 0-gamm)*(1.0 _d 0-gamm)*k1local*k1local |
441 |
& - 4.0 _d 0*k1local*k2local*(1.0 _d 0-2.0 _d 0*gamm) |
442 |
hnew = 0.5 _d 0*( (gamm-1.0 _d 0)*k1local + sqrt(stuff) ) |
443 |
cmick - now determine [CO2*] |
444 |
co2s = diclocal/ |
445 |
& (1.0 _d 0 + (k1local/hnew) + (k1local*k2local/(hnew*hnew))) |
446 |
cmick - return update pH to main routine |
447 |
phlocal = -log10(hnew) |
448 |
|
449 |
c NOW EVALUATE CO32-, carbonate ion concentration |
450 |
c used in determination of calcite compensation depth |
451 |
c Karsten Friis & Mick - Sep 2004 |
452 |
co3local = k1local*k2local*diclocal / |
453 |
& (hnew*hnew + k1local*hnew + k1local*k2local) |
454 |
|
455 |
c --------------------------------------------------------------- |
456 |
c surface pCO2 (following Dickson and Goyet, DOE...) |
457 |
fco2 = co2s/k0local |
458 |
pco2surfloc = fco2/fugflocal |
459 |
|
460 |
C ---------------------------------------------------------------- |
461 |
c Reconvert from mol/kg -> mol/m^3 |
462 |
pt=pt/permil |
463 |
sit=sit/permil |
464 |
ta=ta/permil |
465 |
diclocal=diclocal/permil |
466 |
|
467 |
RETURN |
468 |
END |
469 |
|
470 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
471 |
CBOP |
472 |
subroutine CARBON_COEFFS( |
473 |
I ttemp,stemp, |
474 |
I bi,bj,iMin,iMax,jMin,jMax,myThid) |
475 |
|
476 |
C !DESCRIPTION: |
477 |
C *==========================================================* |
478 |
C | SUBROUTINE CARBON_COEFFS | |
479 |
C | determine coefficients for surface carbon chemistry | |
480 |
C | adapted from OCMIP2: SUBROUTINE CO2CALC | |
481 |
C | mick follows, oct 1999 | |
482 |
C | minor changes to tidy, swd aug 2002 | |
483 |
C *==========================================================* |
484 |
C INPUT |
485 |
C diclocal = total inorganic carbon (mol/m^3) |
486 |
C where 1 T = 1 metric ton = 1000 kg |
487 |
C ta = total alkalinity (eq/m^3) |
488 |
C pt = inorganic phosphate (mol/^3) |
489 |
C sit = inorganic silicate (mol/^3) |
490 |
C t = temperature (degrees C) |
491 |
C s = salinity (PSU) |
492 |
C OUTPUT |
493 |
C IMPORTANT: Some words about units - (JCO, 4/4/1999) |
494 |
C - Models carry tracers in mol/m^3 (on a per volume basis) |
495 |
C - Conversely, this routine, which was written by observationalists |
496 |
C (C. Sabine and R. Key), passes input arguments in umol/kg |
497 |
C (i.e., on a per mass basis) |
498 |
C - I have changed things slightly so that input arguments are in mol/m^3, |
499 |
C - Thus, all input concentrations (diclocal, ta, pt, and st) should be |
500 |
C given in mol/m^3; output arguments "co2star" and "dco2star" |
501 |
C are likewise be in mol/m^3. |
502 |
C |
503 |
C Apr 2011: fix vapour bug (following Bennington) |
504 |
C-------------------------------------------------------------------------- |
505 |
|
506 |
implicit none |
507 |
|
508 |
C === Global variables === |
509 |
#include "SIZE.h" |
510 |
#include "DYNVARS.h" |
511 |
#include "EEPARAMS.h" |
512 |
#include "PARAMS.h" |
513 |
#include "GRID.h" |
514 |
#include "FFIELDS.h" |
515 |
#include "BLING_VARS.h" |
516 |
C == Routine arguments == |
517 |
C ttemp and stemp are local theta and salt arrays |
518 |
C dont really need to pass T and S in, could use theta, salt in |
519 |
C common block in DYNVARS.h, but this way keeps subroutine more |
520 |
C general |
521 |
_RL ttemp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
522 |
_RL stemp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
523 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
524 |
INTEGER myThid |
525 |
CEOP |
526 |
|
527 |
C LOCAL VARIABLES |
528 |
_RL t |
529 |
_RL s |
530 |
_RL tk |
531 |
_RL tk100 |
532 |
_RL tk1002 |
533 |
_RL dlogtk |
534 |
_RL sqrtis |
535 |
_RL sqrts |
536 |
_RL s15 |
537 |
_RL scl |
538 |
_RL s2 |
539 |
_RL invtk |
540 |
_RL is |
541 |
_RL is2 |
542 |
c add Bennington |
543 |
_RL P1atm |
544 |
_RL Rgas |
545 |
_RL RT |
546 |
_RL delta |
547 |
_RL B1 |
548 |
_RL B |
549 |
INTEGER i |
550 |
INTEGER j |
551 |
|
552 |
C..................................................................... |
553 |
C OCMIP note: |
554 |
C Calculate all constants needed to convert between various measured |
555 |
C carbon species. References for each equation are noted in the code. |
556 |
C Once calculated, the constants are |
557 |
C stored and passed in the common block "const". The original version |
558 |
C of this code was based on the code by dickson in Version 2 of |
559 |
C Handbook of Methods C for the Analysis of the Various Parameters of |
560 |
C the Carbon Dioxide System in Seawater , DOE, 1994 (SOP No. 3, p25-26). |
561 |
C.................................................................... |
562 |
|
563 |
do i=imin,imax |
564 |
do j=jmin,jmax |
565 |
if (hFacC(i,j,1,bi,bj).gt.0. _d 0) then |
566 |
t = ttemp(i,j) |
567 |
s = stemp(i,j) |
568 |
C terms used more than once |
569 |
tk = 273.15 _d 0 + t |
570 |
tk100 = tk/100. _d 0 |
571 |
tk1002=tk100*tk100 |
572 |
invtk=1.0 _d 0/tk |
573 |
dlogtk=log(tk) |
574 |
is=19.924 _d 0*s/(1000. _d 0-1.005 _d 0*s) |
575 |
is2=is*is |
576 |
sqrtis=sqrt(is) |
577 |
s2=s*s |
578 |
sqrts=sqrt(s) |
579 |
s15=s**1.5 _d 0 |
580 |
scl=s/1.80655 _d 0 |
581 |
C ----------------------------------------------------------------------- |
582 |
C added by Val Bennington Nov 2010 |
583 |
C Fugacity Factor needed for non-ideality in ocean |
584 |
C ff used for atmospheric correction for water vapor and pressure |
585 |
C Weiss (1974) Marine Chemistry |
586 |
P1atm = 1.01325 _d 0 ! bars |
587 |
Rgas = 83.1451 _d 0 ! bar*cm3/(mol*K) |
588 |
RT = Rgas*tk |
589 |
delta = (57.7 _d 0 - 0.118 _d 0*tk) |
590 |
B1 = -1636.75 _d 0 + 12.0408 _d 0*tk - 0.0327957 _d 0*tk*tk |
591 |
B = B1 + 3.16528 _d 0*tk*tk*tk*(0.00001 _d 0) |
592 |
fugf(i,j,bi,bj) = exp( (B+2. _d 0*delta) * P1atm / RT) |
593 |
C------------------------------------------------------------------------ |
594 |
C f = k0(1-pH2O)*correction term for non-ideality |
595 |
C Weiss & Price (1980, Mar. Chem., 8, 347-359; Eq 13 with table 6 values) |
596 |
ff(i,j,bi,bj) = exp(-162.8301 _d 0 + 218.2968 _d 0/tk100 + |
597 |
& 90.9241 _d 0*log(tk100) - 1.47696 _d 0*tk1002 + |
598 |
& s * (.025695 _d 0 - .025225 _d 0*tk100 + |
599 |
& 0.0049867 _d 0*tk1002)) |
600 |
C------------------------------------------------------------------------ |
601 |
C K0 from Weiss 1974 |
602 |
ak0(i,j,bi,bj) = exp(93.4517 _d 0/tk100 - 60.2409 _d 0 + |
603 |
& 23.3585 _d 0 * log(tk100) + |
604 |
& s * (0.023517 _d 0 - 0.023656 _d 0*tk100 + |
605 |
& 0.0047036 _d 0*tk1002)) |
606 |
C------------------------------------------------------------------------ |
607 |
C k1 = [H][HCO3]/[H2CO3] |
608 |
C k2 = [H][CO3]/[HCO3] |
609 |
C Millero p.664 (1995) using Mehrbach et al. data on seawater scale |
610 |
ak1(i,j,bi,bj)=10.**(-1. _d 0*(3670.7 _d 0*invtk - |
611 |
& 62.008 _d 0 + 9.7944 _d 0*dlogtk - |
612 |
& 0.0118 _d 0 * s + 0.000116 _d 0*s2)) |
613 |
ak2(i,j,bi,bj)=10.**(-1. _d 0*(1394.7 _d 0*invtk+ 4.777 _d 0- |
614 |
& 0.0184 _d 0*s + 0.000118 _d 0*s2)) |
615 |
C------------------------------------------------------------------------ |
616 |
C kb = [H][BO2]/[HBO2] |
617 |
C Millero p.669 (1995) using data from dickson (1990) |
618 |
akb(i,j,bi,bj)=exp((-8966.90 _d 0- 2890.53 _d 0*sqrts - |
619 |
& 77.942 _d 0*s + 1.728 _d 0*s15 - 0.0996 _d 0*s2)*invtk + |
620 |
& (148.0248 _d 0 + 137.1942 _d 0*sqrts + 1.62142 _d 0*s) + |
621 |
& (-24.4344 _d 0 - 25.085 _d 0*sqrts - 0.2474 _d 0*s) * |
622 |
& dlogtk + 0.053105 _d 0*sqrts*tk) |
623 |
C------------------------------------------------------------------------ |
624 |
C k1p = [H][H2PO4]/[H3PO4] |
625 |
C DOE(1994) eq 7.2.20 with footnote using data from Millero (1974) |
626 |
ak1p(i,j,bi,bj) = exp(-4576.752 _d 0*invtk + 115.525 _d 0 - |
627 |
& 18.453 _d 0*dlogtk + |
628 |
& (-106.736 _d 0*invtk + 0.69171 _d 0)*sqrts + |
629 |
& (-0.65643 _d 0*invtk - 0.01844 _d 0)*s) |
630 |
C------------------------------------------------------------------------ |
631 |
C k2p = [H][HPO4]/[H2PO4] |
632 |
C DOE(1994) eq 7.2.23 with footnote using data from Millero (1974)) |
633 |
ak2p(i,j,bi,bj) = exp(-8814.715 _d 0*invtk + 172.0883 _d 0 - |
634 |
& 27.927 _d 0*dlogtk + |
635 |
& (-160.340 _d 0*invtk + 1.3566 _d 0) * sqrts + |
636 |
& (0.37335 _d 0*invtk - 0.05778 _d 0) * s) |
637 |
C------------------------------------------------------------------------ |
638 |
C k3p = [H][PO4]/[HPO4] |
639 |
C DOE(1994) eq 7.2.26 with footnote using data from Millero (1974) |
640 |
ak3p(i,j,bi,bj) = exp(-3070.75 _d 0*invtk - 18.141 _d 0 + |
641 |
& (17.27039 _d 0*invtk + 2.81197 _d 0) * |
642 |
& sqrts + (-44.99486 _d 0*invtk - 0.09984 _d 0) * s) |
643 |
C------------------------------------------------------------------------ |
644 |
C ksi = [H][SiO(OH)3]/[Si(OH)4] |
645 |
C Millero p.671 (1995) using data from Yao and Millero (1995) |
646 |
aksi(i,j,bi,bj) = exp(-8904.2 _d 0*invtk + 117.385 _d 0 - |
647 |
& 19.334 _d 0*dlogtk + |
648 |
& (-458.79 _d 0*invtk + 3.5913 _d 0) * sqrtis + |
649 |
& (188.74 _d 0*invtk - 1.5998 _d 0) * is + |
650 |
& (-12.1652 _d 0*invtk + 0.07871 _d 0) * is2 + |
651 |
& log(1.0 _d 0-0.001005 _d 0*s)) |
652 |
C------------------------------------------------------------------------ |
653 |
C kw = [H][OH] |
654 |
C Millero p.670 (1995) using composite data |
655 |
akw(i,j,bi,bj) = exp(-13847.26 _d 0*invtk + 148.9652 _d 0 - |
656 |
& 23.6521 _d 0*dlogtk + |
657 |
& (118.67 _d 0*invtk - 5.977 _d 0 + 1.0495 _d 0 * dlogtk) |
658 |
& * sqrts - 0.01615 _d 0 * s) |
659 |
C------------------------------------------------------------------------ |
660 |
C ks = [H][SO4]/[HSO4] |
661 |
C dickson (1990, J. chem. Thermodynamics 22, 113) |
662 |
aks(i,j,bi,bj)=exp(-4276.1 _d 0*invtk + 141.328 _d 0 - |
663 |
& 23.093 _d 0*dlogtk + |
664 |
& (-13856. _d 0*invtk + 324.57 _d 0 - 47.986 _d 0*dlogtk)*sqrtis+ |
665 |
& (35474. _d 0*invtk - 771.54 _d 0 + 114.723 _d 0*dlogtk)*is - |
666 |
& 2698. _d 0*invtk*is**1.5 _d 0 + 1776. _d 0*invtk*is2 + |
667 |
& log(1.0 _d 0 - 0.001005 _d 0*s)) |
668 |
C------------------------------------------------------------------------ |
669 |
C kf = [H][F]/[HF] |
670 |
C dickson and Riley (1979) -- change pH scale to total |
671 |
akf(i,j,bi,bj)=exp(1590.2 _d 0*invtk - 12.641 _d 0 + |
672 |
& 1.525 _d 0*sqrtis + log(1.0 _d 0 - 0.001005 _d 0*s) + |
673 |
& log(1.0 _d 0 + (0.1400 _d 0/96.062 _d 0)*(scl)/aks(i,j,bi,bj))) |
674 |
C------------------------------------------------------------------------ |
675 |
C Calculate concentrations for borate, sulfate, and fluoride |
676 |
C Uppstrom (1974) |
677 |
bt(i,j,bi,bj) = 0.000232 _d 0 * scl/10.811 _d 0 |
678 |
C Morris & Riley (1966) |
679 |
st(i,j,bi,bj) = 0.14 _d 0 * scl/96.062 _d 0 |
680 |
C Riley (1965) |
681 |
ft(i,j,bi,bj) = 0.000067 _d 0 * scl/18.9984 _d 0 |
682 |
C------------------------------------------------------------------------ |
683 |
else |
684 |
c add Bennington |
685 |
fugf(i,j,bi,bj)=0. _d 0 |
686 |
ff(i,j,bi,bj)=0. _d 0 |
687 |
ak0(i,j,bi,bj)= 0. _d 0 |
688 |
ak1(i,j,bi,bj)= 0. _d 0 |
689 |
ak2(i,j,bi,bj)= 0. _d 0 |
690 |
akb(i,j,bi,bj)= 0. _d 0 |
691 |
ak1p(i,j,bi,bj) = 0. _d 0 |
692 |
ak2p(i,j,bi,bj) = 0. _d 0 |
693 |
ak3p(i,j,bi,bj) = 0. _d 0 |
694 |
aksi(i,j,bi,bj) = 0. _d 0 |
695 |
akw(i,j,bi,bj) = 0. _d 0 |
696 |
aks(i,j,bi,bj)= 0. _d 0 |
697 |
akf(i,j,bi,bj)= 0. _d 0 |
698 |
bt(i,j,bi,bj) = 0. _d 0 |
699 |
st(i,j,bi,bj) = 0. _d 0 |
700 |
ft(i,j,bi,bj) = 0. _d 0 |
701 |
endif |
702 |
end do |
703 |
end do |
704 |
|
705 |
RETURN |
706 |
END |
707 |
|
708 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
709 |
|
710 |
CBOP |
711 |
subroutine CARBON_COEFFS_PRESSURE_DEP( |
712 |
I ttemp,stemp, |
713 |
I bi,bj,iMin,iMax,jMin,jMax, |
714 |
I Klevel,myThid) |
715 |
|
716 |
C !DESCRIPTION: |
717 |
C *==========================================================* |
718 |
C | SUBROUTINE CARBON_COEFFS | |
719 |
C | determine coefficients for surface carbon chemistry | |
720 |
C | adapted from OCMIP2: SUBROUTINE CO2CALC | |
721 |
C | mick follows, oct 1999 | |
722 |
C | minor changes to tidy, swd aug 2002 | |
723 |
C | MODIFIED FOR PRESSURE DEPENDENCE | |
724 |
C | Karsten Friis and Mick Follows 2004 | |
725 |
C *==========================================================* |
726 |
C INPUT |
727 |
C diclocal = total inorganic carbon (mol/m^3) |
728 |
C where 1 T = 1 metric ton = 1000 kg |
729 |
C ta = total alkalinity (eq/m^3) |
730 |
C pt = inorganic phosphate (mol/^3) |
731 |
C sit = inorganic silicate (mol/^3) |
732 |
C t = temperature (degrees C) |
733 |
C s = salinity (PSU) |
734 |
C OUTPUT |
735 |
C IMPORTANT: Some words about units - (JCO, 4/4/1999) |
736 |
C - Models carry tracers in mol/m^3 (on a per volume basis) |
737 |
C - Conversely, this routine, which was written by observationalists |
738 |
C (C. Sabine and R. Key), passes input arguments in umol/kg |
739 |
C (i.e., on a per mass basis) |
740 |
C - I have changed things slightly so that input arguments are in mol/m^3, |
741 |
C - Thus, all input concentrations (diclocal, ta, pt, and st) should be |
742 |
C given in mol/m^3; output arguments "co2star" and "dco2star" |
743 |
C are likewise be in mol/m^3. |
744 |
C |
745 |
C NOW INCLUDES: |
746 |
C PRESSURE DEPENDENCE of K1, K2, solubility product of calcite |
747 |
C based on Takahashi, GEOSECS Atlantic Report, Vol. 1 (1981) |
748 |
C |
749 |
C-------------------------------------------------------------------------- |
750 |
|
751 |
implicit none |
752 |
|
753 |
C === Global variables === |
754 |
#include "SIZE.h" |
755 |
#include "DYNVARS.h" |
756 |
#include "EEPARAMS.h" |
757 |
#include "PARAMS.h" |
758 |
#include "GRID.h" |
759 |
#include "FFIELDS.h" |
760 |
#include "BLING_VARS.h" |
761 |
C == Routine arguments == |
762 |
C ttemp and stemp are local theta and salt arrays |
763 |
C dont really need to pass T and S in, could use theta, salt in |
764 |
C common block in DYNVARS.h, but this way keeps subroutine more |
765 |
C general |
766 |
_RL ttemp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
767 |
_RL stemp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
768 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
769 |
C K is depth index |
770 |
INTEGER Klevel |
771 |
INTEGER myThid |
772 |
CEOP |
773 |
|
774 |
C LOCAL VARIABLES |
775 |
_RL t |
776 |
_RL s |
777 |
_RL tk |
778 |
_RL tk100 |
779 |
_RL tk1002 |
780 |
_RL dlogtk |
781 |
_RL sqrtis |
782 |
_RL sqrts |
783 |
_RL s15 |
784 |
_RL scl |
785 |
_RL s2 |
786 |
_RL invtk |
787 |
_RL is |
788 |
_RL is2 |
789 |
c add Bennington |
790 |
_RL P1atm |
791 |
_RL Rgas |
792 |
_RL RT |
793 |
_RL delta |
794 |
_RL B1 |
795 |
_RL B |
796 |
INTEGER i |
797 |
INTEGER j |
798 |
INTEGER k |
799 |
_RL bdepth |
800 |
_RL cdepth |
801 |
_RL pressc |
802 |
_RL Ksp_T_Calc |
803 |
_RL Ksp_T_Arag |
804 |
_RL xvalue |
805 |
_RL zdum |
806 |
_RL tmpa1 |
807 |
_RL tmpa2 |
808 |
_RL tmpa3 |
809 |
_RL logKspc |
810 |
_RL logKspa |
811 |
_RL dv |
812 |
_RL dk |
813 |
_RL pfactor |
814 |
_RL bigR |
815 |
|
816 |
C..................................................................... |
817 |
C OCMIP note: |
818 |
C Calculate all constants needed to convert between various measured |
819 |
C carbon species. References for each equation are noted in the code. |
820 |
C Once calculated, the constants are |
821 |
C stored and passed in the common block "const". The original version |
822 |
C of this code was based on the code by dickson in Version 2 of |
823 |
C Handbook of Methods C for the Analysis of the Various Parameters of |
824 |
C the Carbon Dioxide System in Seawater , DOE, 1994 (SOP No. 3, p25-26). |
825 |
C.................................................................... |
826 |
|
827 |
c determine pressure (bar) from depth |
828 |
c 1 BAR at z=0m (atmos pressure) |
829 |
c use UPPER surface of cell so top layer pressure = 0 bar |
830 |
c for surface exchange coeffs |
831 |
|
832 |
cmick.............................. |
833 |
c write(6,*)'Klevel ',klevel |
834 |
|
835 |
bdepth = 0.0d0 |
836 |
cdepth = 0.0d0 |
837 |
pressc = 1.0d0 |
838 |
do k = 1,Klevel |
839 |
cdepth = bdepth + 0.5d0*drF(k) |
840 |
bdepth = bdepth + drF(k) |
841 |
pressc = 1.0d0 + 0.1d0*cdepth |
842 |
end do |
843 |
|
844 |
do i=imin,imax |
845 |
do j=jmin,jmax |
846 |
if (hFacC(i,j,Klevel,bi,bj).gt.0.d0) then |
847 |
t = ttemp(i,j) |
848 |
s = stemp(i,j) |
849 |
C terms used more than once |
850 |
tk = 273.15 + t |
851 |
tk100 = tk/100.0 |
852 |
tk1002=tk100*tk100 |
853 |
invtk=1.0/tk |
854 |
dlogtk=log(tk) |
855 |
is=19.924*s/(1000.-1.005*s) |
856 |
is2=is*is |
857 |
sqrtis=sqrt(is) |
858 |
s2=s*s |
859 |
sqrts=sqrt(s) |
860 |
s15=s**1.5 |
861 |
scl=s/1.80655 |
862 |
|
863 |
C----------------------------------------------------------------------- |
864 |
C Ariane adding this - copied from non-pressure-dep case |
865 |
C added by Val Bennington Nov 2010 |
866 |
C Fugacity Factor needed for non-ideality in ocean |
867 |
C ff used for atmospheric correction for water vapor and pressure |
868 |
C Weiss (1974) Marine Chemistry |
869 |
P1atm = 1.01325 _d 0 ! bars |
870 |
Rgas = 83.1451 _d 0 ! bar*cm3/(mol*K) |
871 |
RT = Rgas*tk |
872 |
delta = (57.7 _d 0 - 0.118 _d 0*tk) |
873 |
B1 = -1636.75 _d 0 + 12.0408 _d 0*tk - 0.0327957 _d 0*tk*tk |
874 |
B = B1 + 3.16528 _d 0*tk*tk*tk*(0.00001 _d 0) |
875 |
fugf(i,j,bi,bj) = exp( (B+2. _d 0*delta) * P1atm / RT) |
876 |
|
877 |
C------------------------------------------------------------------------ |
878 |
C f = k0(1-pH2O)*correction term for non-ideality |
879 |
C Weiss & Price (1980, Mar. Chem., 8, 347-359; Eq 13 with table 6 values) |
880 |
ff(i,j,bi,bj) = exp(-162.8301 + 218.2968/tk100 + |
881 |
& 90.9241*log(tk100) - 1.47696*tk1002 + |
882 |
& s * (.025695 - .025225*tk100 + |
883 |
& 0.0049867*tk1002)) |
884 |
C------------------------------------------------------------------------ |
885 |
C K0 from Weiss 1974 |
886 |
ak0(i,j,bi,bj) = exp(93.4517/tk100 - 60.2409 + |
887 |
& 23.3585 * log(tk100) + |
888 |
& s * (0.023517 - 0.023656*tk100 + |
889 |
& 0.0047036*tk1002)) |
890 |
C------------------------------------------------------------------------ |
891 |
C k1 = [H][HCO3]/[H2CO3] |
892 |
C k2 = [H][CO3]/[HCO3] |
893 |
C Millero p.664 (1995) using Mehrbach et al. data on seawater scale |
894 |
ak1(i,j,bi,bj)=10**(-1*(3670.7*invtk - |
895 |
& 62.008 + 9.7944*dlogtk - |
896 |
& 0.0118 * s + 0.000116*s2)) |
897 |
ak2(i,j,bi,bj)=10**(-1*(1394.7*invtk + 4.777 - |
898 |
& 0.0184*s + 0.000118*s2)) |
899 |
C NOW PRESSURE DEPENDENCE: |
900 |
c Following Takahashi (1981) GEOSECS report - quoting Culberson and |
901 |
c Pytkowicz (1968) |
902 |
c pressc = pressure in bars |
903 |
ak1(i,j,bi,bj) = ak1(i,j,bi,bj)* |
904 |
& exp( (24.2-0.085*t)*(pressc-1.0)/(83.143*tk) ) |
905 |
c FIRST GO FOR K2: According to GEOSECS (1982) report |
906 |
c ak2(i,j,bi,bj) = ak2(i,j,bi,bj)* |
907 |
c & exp( (26.4-0.040*t)*(pressc-1.0)/(83.143*tk) ) |
908 |
c SECOND GO FOR K2: corrected coeff according to CO2sys documentation |
909 |
c E. Lewis and D. Wallace (1998) ORNL/CDIAC-105 |
910 |
ak2(i,j,bi,bj) = ak2(i,j,bi,bj)* |
911 |
& exp( (16.4-0.040*t)*(pressc-1.0)/(83.143*tk) ) |
912 |
C------------------------------------------------------------------------ |
913 |
C kb = [H][BO2]/[HBO2] |
914 |
C Millero p.669 (1995) using data from dickson (1990) |
915 |
akb(i,j,bi,bj)=exp((-8966.90 - 2890.53*sqrts - 77.942*s + |
916 |
& 1.728*s15 - 0.0996*s2)*invtk + |
917 |
& (148.0248 + 137.1942*sqrts + 1.62142*s) + |
918 |
& (-24.4344 - 25.085*sqrts - 0.2474*s) * |
919 |
& dlogtk + 0.053105*sqrts*tk) |
920 |
C Mick and Karsten - Dec 04 |
921 |
C ADDING pressure dependence based on Millero (1995), p675 |
922 |
C with additional info from CO2sys documentation (E. Lewis and |
923 |
C D. Wallace, 1998 - see endnotes for commentary on Millero, 95) |
924 |
bigR = 83.145 |
925 |
dv = -29.48 + 0.1622*t + 2.608d-3*t*t |
926 |
dk = -2.84d-3 |
927 |
pfactor = - (dv/(bigR*tk))*pressc |
928 |
& + (0.5*dk/(bigR*tk))*pressc*pressc |
929 |
akb(i,j,bi,bj) = akb(i,j,bi,bj)*exp(pfactor) |
930 |
C------------------------------------------------------------------------ |
931 |
C k1p = [H][H2PO4]/[H3PO4] |
932 |
C DOE(1994) eq 7.2.20 with footnote using data from Millero (1974) |
933 |
ak1p(i,j,bi,bj) = exp(-4576.752*invtk + 115.525 - |
934 |
& 18.453*dlogtk + |
935 |
& (-106.736*invtk + 0.69171)*sqrts + |
936 |
& (-0.65643*invtk - 0.01844)*s) |
937 |
C------------------------------------------------------------------------ |
938 |
C k2p = [H][HPO4]/[H2PO4] |
939 |
C DOE(1994) eq 7.2.23 with footnote using data from Millero (1974)) |
940 |
ak2p(i,j,bi,bj) = exp(-8814.715*invtk + 172.0883 - |
941 |
& 27.927*dlogtk + |
942 |
& (-160.340*invtk + 1.3566) * sqrts + |
943 |
& (0.37335*invtk - 0.05778) * s) |
944 |
C------------------------------------------------------------------------ |
945 |
C k3p = [H][PO4]/[HPO4] |
946 |
C DOE(1994) eq 7.2.26 with footnote using data from Millero (1974) |
947 |
ak3p(i,j,bi,bj) = exp(-3070.75*invtk - 18.141 + |
948 |
& (17.27039*invtk + 2.81197) * |
949 |
& sqrts + (-44.99486*invtk - 0.09984) * s) |
950 |
C------------------------------------------------------------------------ |
951 |
C ksi = [H][SiO(OH)3]/[Si(OH)4] |
952 |
C Millero p.671 (1995) using data from Yao and Millero (1995) |
953 |
aksi(i,j,bi,bj) = exp(-8904.2*invtk + 117.385 - |
954 |
& 19.334*dlogtk + |
955 |
& (-458.79*invtk + 3.5913) * sqrtis + |
956 |
& (188.74*invtk - 1.5998) * is + |
957 |
& (-12.1652*invtk + 0.07871) * is2 + |
958 |
& log(1.0-0.001005*s)) |
959 |
C------------------------------------------------------------------------ |
960 |
C kw = [H][OH] |
961 |
C Millero p.670 (1995) using composite data |
962 |
akw(i,j,bi,bj) = exp(-13847.26*invtk + 148.9652 - |
963 |
& 23.6521*dlogtk + |
964 |
& (118.67*invtk - 5.977 + 1.0495 * dlogtk) * |
965 |
& sqrts - 0.01615 * s) |
966 |
C------------------------------------------------------------------------ |
967 |
C ks = [H][SO4]/[HSO4] |
968 |
C dickson (1990, J. chem. Thermodynamics 22, 113) |
969 |
aks(i,j,bi,bj)=exp(-4276.1*invtk + 141.328 - |
970 |
& 23.093*dlogtk + |
971 |
& (-13856*invtk + 324.57 - 47.986*dlogtk)*sqrtis + |
972 |
& (35474*invtk - 771.54 + 114.723*dlogtk)*is - |
973 |
& 2698*invtk*is**1.5 + 1776*invtk*is2 + |
974 |
& log(1.0 - 0.001005*s)) |
975 |
C------------------------------------------------------------------------ |
976 |
C kf = [H][F]/[HF] |
977 |
C dickson and Riley (1979) -- change pH scale to total |
978 |
akf(i,j,bi,bj)=exp(1590.2*invtk - 12.641 + 1.525*sqrtis + |
979 |
& log(1.0 - 0.001005*s) + |
980 |
& log(1.0 + (0.1400/96.062)*(scl)/aks(i,j,bi,bj))) |
981 |
C------------------------------------------------------------------------ |
982 |
C Calculate concentrations for borate, sulfate, and fluoride |
983 |
C Uppstrom (1974) |
984 |
bt(i,j,bi,bj) = 0.000232 * scl/10.811 |
985 |
C Morris & Riley (1966) |
986 |
st(i,j,bi,bj) = 0.14 * scl/96.062 |
987 |
C Riley (1965) |
988 |
ft(i,j,bi,bj) = 0.000067 * scl/18.9984 |
989 |
C------------------------------------------------------------------------ |
990 |
C solubility product for calcite |
991 |
caxx and aragonite |
992 |
C |
993 |
c Following Takahashi (1982) GEOSECS handbook |
994 |
C NOT SURE THIS IS WORKING??? |
995 |
C Ingle et al. (1973) |
996 |
c Ksp_T_Calc = ( -34.452 - 39.866*(s**0.333333) |
997 |
c & + 110.21*log(s) - 7.5752d-6 * (tk**2.0) |
998 |
c & ) * 1.0d-7 |
999 |
c with pressure dependence Culberson and Pytkowicz (1968) |
1000 |
c xvalue = (36-0.20*t)*(pressc-1.0)/(83.143*tk) |
1001 |
c Ksp_TP_Calc(i,j,bi,bj) = Ksp_T_Calc*exp(xvalue) |
1002 |
c |
1003 |
c |
1004 |
C Following Mucci (1983) - from Zeebe/Wolf-Gladrow equic.m |
1005 |
tmpa1 = - 171.9065 - (0.077993*tk) + (2839.319/tk) |
1006 |
& + (71.595*log10(tk)) |
1007 |
tmpa2 = +(-0.77712 + (0.0028426*tk) + (178.34/tk) )*sqrts |
1008 |
tmpa3 = -(0.07711*s) + (0.0041249*s15) |
1009 |
logKspc = tmpa1 + tmpa2 + tmpa3 |
1010 |
Ksp_T_Calc = 10.0**logKspc |
1011 |
C |
1012 |
tmpa1 = - 171.945 - (0.077993*tk) + (2903.293/tk) |
1013 |
& + (71.595*log10(tk)) |
1014 |
tmpa2 = +(-0.068393 + (0.0017276*tk) + (88.135/tk) )*sqrts |
1015 |
tmpa3 = -(0.10018*s) + (0.0059415*s15) |
1016 |
logKspa = tmpa1 + tmpa2 + tmpa3 |
1017 |
Ksp_T_Arag = 10.0**logKspa |
1018 |
|
1019 |
c with pressure dependence Culberson and Pytkowicz (1968) |
1020 |
c xvalue = (36.0-0.20*t)*(pressc-1.0)/(83.143*tk) |
1021 |
c Ksp_TP_Calc(i,j,bi,bj) = Ksp_T_Calc*exp(xvalue) |
1022 |
|
1023 |
c alternative pressure depdendence |
1024 |
c following Millero (1995) but using info from Appendix A11 of |
1025 |
c Zeebe and Wolf-Gladrow (2001) book |
1026 |
c dv = -48.6 - 0.5304*t |
1027 |
c dk = -11.76d-3 - 0.3692*t |
1028 |
c xvalue = - (dv/(bigR*tk))*pressc |
1029 |
c & + (0.5*dk/(bigR*tk))*pressc*pressc |
1030 |
c Ksp_TP_Calc(i,j,bi,bj) = Ksp_T_Calc*exp(xvalue) |
1031 |
|
1032 |
c alternative pressure dependence from Ingle (1975) |
1033 |
|
1034 |
zdum = (pressc*10.0d0 - 10.0d0)/10.0d0 |
1035 |
xvalue = ( (48.8d0 - 0.53d0*t)*zdum |
1036 |
& + (-0.00588d0 + 0.0001845d0*t)*zdum*zdum) |
1037 |
& / (188.93d0*(t + 273.15d0)) |
1038 |
|
1039 |
Ksp_TP_Calc(i,j,bi,bj) = Ksp_T_Calc*10**(xvalue) |
1040 |
Ksp_TP_Arag(i,j,bi,bj) = Ksp_T_Arag*10**(xvalue) |
1041 |
|
1042 |
|
1043 |
C------------------------------------------------------------------------ |
1044 |
else |
1045 |
fugf(i,j,bi,bj)=0.d0 |
1046 |
ff(i,j,bi,bj)=0.d0 |
1047 |
ak0(i,j,bi,bj)= 0.d0 |
1048 |
ak1(i,j,bi,bj)= 0.d0 |
1049 |
ak2(i,j,bi,bj)= 0.d0 |
1050 |
akb(i,j,bi,bj)= 0.d0 |
1051 |
ak1p(i,j,bi,bj) = 0.d0 |
1052 |
ak2p(i,j,bi,bj) = 0.d0 |
1053 |
ak3p(i,j,bi,bj) = 0.d0 |
1054 |
aksi(i,j,bi,bj) = 0.d0 |
1055 |
akw(i,j,bi,bj) = 0.d0 |
1056 |
aks(i,j,bi,bj)= 0.d0 |
1057 |
akf(i,j,bi,bj)= 0.d0 |
1058 |
bt(i,j,bi,bj) = 0.d0 |
1059 |
st(i,j,bi,bj) = 0.d0 |
1060 |
ft(i,j,bi,bj) = 0.d0 |
1061 |
Ksp_TP_Calc(i,j,bi,bj) = 0.d0 |
1062 |
Ksp_TP_Arag(i,j,bi,bj) = 0.d0 |
1063 |
endif |
1064 |
end do |
1065 |
end do |
1066 |
|
1067 |
return |
1068 |
end |