/[MITgcm]/MITgcm/pkg/dic/carbon_chem.F
ViewVC logotype

Contents of /MITgcm/pkg/dic/carbon_chem.F

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


Revision 1.13 - (show annotations) (download)
Wed Apr 9 22:13:15 2008 UTC (16 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint62, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.12: +14 -7 lines
add argument "myThid" where missing

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

  ViewVC Help
Powered by ViewVC 1.1.22