/[MITgcm]/MITgcm_contrib/darwin2/pkg/darwin/dic_carbon_chem.F
ViewVC logotype

Annotation of /MITgcm_contrib/darwin2/pkg/darwin/dic_carbon_chem.F

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


Revision 1.3 - (hide annotations) (download)
Mon May 9 20:17:39 2011 UTC (14 years, 2 months ago) by stephd
Branch: MAIN
CVS Tags: ctrb_darwin2_ckpt64k_20130723, ctrb_darwin2_ckpt63l_20120405, ctrb_darwin2_ckpt64h_20130528, ctrb_darwin2_ckpt64m_20130820, ctrb_darwin2_ckpt64f_20130405, ctrb_darwin2_ckpt63f_20111201, ctrb_darwin2_ckpt64a_20121116, ctrb_darwin2_ckpt64n_20130826, ctrb_darwin2_ckpt62y_20110526, ctrb_darwin2_ckpt64i_20130622, ctrb_darwin2_ckpt62x_20110513, ctrb_darwin2_ckpt63o_20120629, ctrb_darwin2_ckpt64e_20130305, ctrb_darwin2_ckpt63c_20111011, ctrb_darwin2_ckpt63i_20120124, ctrb_darwin2_ckpt63m_20120506, ctrb_darwin2_ckpt63s_20120908, ctrb_darwin2_ckpt63e_20111107, ctrb_darwin2_ckpt63b_20110830, ctrb_darwin2_ckpt63j_20120217, ctrb_darwin2_ckpt63r_20120817, ctrb_darwin2_ckpt64g_20130503, ctrb_darwin2_ckpt64l_20130806, ctrb_darwin2_ckpt63g_20111220, ctrb_darwin2_ckpt64c_20130120, ctrb_darwin2_ckpt63a_20110804, ctrb_darwin2_ckpt64j_20130704, ctrb_darwin2_ckpt63h_20111230, ctrb_darwin2_ckpt63p_20120707, ctrb_darwin2_ckpt63d_20111107, ctrb_darwin2_ckpt63q_20120731, ctrb_darwin2_ckpt63_20110728, ctrb_darwin2_ckpt64b_20121224, ctrb_darwin2_ckpt64d_20130219, ctrb_darwin2_ckpt64_20121012, ctrb_darwin2_ckpt63n_20120604, ctrb_darwin2_ckpt63k_20120317, ctrb_darwin2_ckpt62z_20110622
Changes since 1.2: +4 -4 lines
o add extra bounds to dic/temp/salt so pH solver doesn't blow up

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

  ViewVC Help
Powered by ViewVC 1.1.22