/[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.1 - (hide annotations) (download)
Wed Apr 13 18:56:24 2011 UTC (14 years, 3 months ago) by jahn
Branch: MAIN
CVS Tags: ctrb_darwin2_ckpt62v_20110413, ctrb_darwin2_baseline
darwin2 initial checkin

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

  ViewVC Help
Powered by ViewVC 1.1.22