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