-
Notifications
You must be signed in to change notification settings - Fork 6
/
Geochemistry.jl
2620 lines (2276 loc) · 121 KB
/
Geochemistry.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
## --- Calculate Eu*
"""
```julia
eustar(Nd::Number, Sm::Number, Gd::Number, Tb::Number)
```
Calculate expected europium concentration, Eu*, based on abundance of
adjacent rare earths.
Full four-element log-linear interpolation, assuming 3+ ionic radii and the
chondritic abundances of Sun and McDonough 1989 (doi: 10.1144/gsl.sp.1989.042.01.19)
"""
function eustar(Nd::Number, Sm::Number, Gd::Number, Tb::Number)
# Ionic radii, in pm [Tb, Gd, Sm, Nd]
r = [106.3, 107.8, 109.8, 112.3] # or x = [1, 2, 4, 6]
# Normalize to chondrite
y = log.([Tb/0.0374, Gd/0.2055, Sm/0.1530, Nd/0.4670])
notnan = .!isnan.(y)
# Make sure we're interpolating and not extrapolating
if any(view(notnan, 1:2)) && any(view(notnan, 3:4))
# Fit a straight line through the chondrite-normalized values
x = r[notnan]
(a,b) = hcat(fill!(similar(x), 1), x) \ y[notnan]
# De-dormalize output for Eu, interpolating at r = 108.7 pm or x = 3
eu_interp = 0.0580*exp(a + b*108.7)
else
eu_interp = NaN
end
return eu_interp
end
"""
```julia
eustar(Sm::Number, Gd::Number)
```
Calculate expected europium concentration, Eu*, based on abundance of
adjacent rare earths.
Simple geometric mean interpolation from Sm and Gd alone, assuming the chondritic
abundances of Sun and McDonough 1989 (doi: 10.1144/gsl.sp.1989.042.01.19), that is
Eu* = `0.0580*sqrt(Sm/0.1530 * Gd/0.2055)`
"""
function eustar(Sm::Number, Gd::Number)
# Geometric mean in regular space is equal to the arithmetic mean in log space. Fancy that!
return 0.0580*sqrt(Sm/0.1530 * Gd/0.2055)
end
export eustar
## --- CIPW norm
"""
```julia
cipw_norm(SiO2, TiO2, Al2O3, Fe2O3, FeO, MnO, MgO, CaO, Na2O, K2O, P2O5)
```
Returns
```
quartz, orthoclase, plagioclase, corundum, nepheline, diopside, orthopyroxene, olivine, magnetite, ilmenite, apatite
```
"""
function cipw_norm(SiO2, TiO2, Al2O3, Fe2O3, FeO, MnO, MgO, CaO, Na2O, K2O, P2O5)
SiO2 /= 60.0843
TiO2 /= 79.8988
Al2O3 /= 101.9613
Fe2O3 /= 159.6922
FeO /= 71.8464
MnO /= 70.9374
MgO /= 40.3044
CaO /= 56.0794
Na2O /= 61.9789
K2O /= 94.1960
P2O5 /= 141.9445
FeO = nanadd(FeO, MnO)
CaO -= 3.333333333333333 * P2O5
apatite = 0.6666666666666666 * P2O5
# P2O5 = 0
FeO -= TiO2
ilmenite = TiO2
FeO -= Fe2O3
magnetite = Fe2O3
# Fe2O3 = 0
Al2O3 -= K2O
orthoclase = K2O
# K2O = 0
Al2O3 -= Na2O
albite = Na2O
if CaO > Al2O3
CaO -= Al2O3
anorthite = Al2O3
Al2O3 = 0
else
Al2O3 -= CaO
anorthite = CaO
CaO = 0
end
if Al2O3 > 0
corundum = Al2O3
Al2O3 = 0
else
corundum = 0
end
Mg′ = MgO / (MgO + FeO)
FMO = FeO + MgO
FMO_weight = (Mg′*40.3044)+((1-Mg′)*71.8464)
if CaO > 0
FMO -= CaO
diopside = CaO
else
diopside = 0
end
orthopyroxene = FMO
pSi1 = 6orthoclase + 6albite + 2anorthite + 2diopside + orthopyroxene
if pSi1 < SiO2
quartz = SiO2 - pSi1
nepheline = 0
olivine = 0
else
quartz = 0
pSi2 = 6orthoclase + 6albite + 2anorthite + 2diopside
pSi3 = SiO2 - pSi2
if FMO > 2pSi3
orthopyroxene = 0
olivine = FMO
FMO = 0
pSi4 = 6orthoclase + 2anorthite + 2diopside + 0.5olivine
pSi5 = SiO2 - pSi4
Albite = (pSi5-(2*Na2O))/4
nepheline = Na2O-Albite
else
nepheline = 0
orthopyroxene = 2pSi3 - FMO
olivine = FMO - pSi3
end
end
orthoclase *= 2
nepheline *= 2
albite *= 2
An′ = anorthite/(anorthite+albite)
plag_weight = (An′*278.2093)+((1-An′)*262.2230)
plagioclase = albite+anorthite
quartz *= 60.0843
orthoclase *= 278.3315
plagioclase *= plag_weight
corundum *= 101.9613
nepheline *= 142.0544
diopside *= (172.248 + FMO_weight)
orthopyroxene *= (60.0843 + FMO_weight)
olivine *= (60.0843 + 2FMO_weight)
magnetite *= 231.5386
ilmenite *= 151.7452
apatite *= 504.3152
return (quartz=quartz, orthoclase=orthoclase, plagioclase=plagioclase,
corundum=corundum, nepheline=nepheline, diopside=diopside,
orthopyroxene=orthopyroxene, olivine=olivine, magnetite=magnetite,
ilmenite=ilmenite, apatite=apatite)
end
# export cipw_norm
## --- Fe oxide conversions
"""
```julia
feoconversion(FeO::Number=NaN, Fe2O3::Number=NaN, FeOT::Number=NaN, Fe2O3T::Number=NaN)
```
Compiles data from FeO, Fe2O3, FeOT, and Fe2O3T into a single FeOT value.
"""
function feoconversion(FeO::Number=NaN, Fe2O3::Number=NaN, FeOT::Number=NaN, Fe2O3T::Number=NaN)
# To convert from Fe2O3 wt % to FeO wt %, multiply by
conversionfactor = (55.845+15.999) / (55.845+1.5*15.999)
# If FeOT or Fe2O3T already exists, use that
if isnan(FeOT)
if isnan(Fe2O3T)
if isnan(Fe2O3)
FeOT = FeO
elseif isnan(FeO)
FeOT = Fe2O3*conversionfactor
else
FeOT = Fe2O3*conversionfactor + FeO
end
else
FeOT=Fe2O3T*conversionfactor
end
end
return FeOT
end
export feoconversion
## --- Oxide conversions
function fillifnan!(dest::AbstractArray, source::AbstractArray)
@inbounds for i in eachindex(dest, source)
if isnan(dest[i]) && !isnan(source[i])
dest[i] = source[i]
end
end
return dest
end
function fillifnan!(dest::AbstractArray, source::AbstractArray, factor::Number)
@inbounds for i in eachindex(dest, source)
if isnan(dest[i]) && !isnan(source[i])
dest[i] = source[i] * factor
end
end
return dest
end
function nannegative!(a::AbstractArray)
@inbounds for i in eachindex(a)
if a[i] < 0
a[i] = NaN
end
end
return a
end
"""
```julia
converted_dataset = oxideconversion(dataset::Union{Dict,NamedTuple}; unitratio::Number=10000)
```
As `oxideconversion!`, but returning a copy rather than modifying in-place
"""
oxideconversion(ds::Union{Dict,NamedTuple}; kwargs...) = oxideconversion!(deepcopy(ds); kwargs...)
export oxideconversion
"""
```julia
oxideconversion!(dataset::Dict; unitratio::Number=10000)
```
Convert major elements (Ti, Al, etc.) into corresponding oxides (TiO2, Al2O3, ...) in place if extant.
If metals are expected as PPM, set unitratio=10000 (default); if metals are as wt%,
set unitratio = 1
See also `oxideconversion`, c.f. `metalconversion!`
"""
function oxideconversion!(dataset::NamedTuple; unitratio::Number=10000)
# List of elements to convert
source = (:Si, :Ti, :Al, :Fe, :Fe, :Mg, :Ca, :Mn, :Li, :Na, :K, :P, :Cr, :Ni, :Co, :C, :S, :H)
dest = (:SiO2, :TiO2, :Al2O3, :FeOT, :Fe2O3T, :MgO, :CaO, :MnO, :Li2O, :Na2O, :K2O, :P2O5, :Cr2O3, :NiO, :CoO, :CO2, :SO3, :H2O)
conversionfactor = (2.13932704290547,1.66847584248889,1.88944149488507,1.28648836426407,1.42973254639611,1.65825961736268,1.39919258253823,1.29121895771597,2.1526657060518732,1.34795912485574,1.20459963614796,2.29133490474735,1.46154369861159,1.27258582901258,1.27147688434143,3.66405794688203,2.4970991890205863,8.93601190476191)
# If source field exists, fill in destination from source
for i ∈ eachindex(source)
if haskey(dataset, source[i])
if haskey(dataset, dest[i]) # If destination field doesn't exist, make it.
oxide, metal = dataset[dest[i]], dataset[source[i]]
fillifnan!(oxide, metal, conversionfactor[i]/unitratio)
end
end
end
return dataset
end
oxideconversion!(ds::Dict; kwargs...) = (oxideconversion!(TupleDataset(ds); kwargs...); ds)
export oxideconversion!
"""
```julia
converted_dataset = metalconversion(dataset::Union{Dict,NamedTuple}; unitratio::Number=10000)
```
As `metalconversion!`, but returning a copy rather than modifying in-place
"""
metalconversion(ds::Union{Dict,NamedTuple}; kwargs...) = metalconversion!(copy(ds); kwargs...)
export metalconversion
"""
```julia
dataset = metalconversion!(dataset::Union{Dict,NamedTuple}; unitratio::Number=10000)
```
Convert minor element oxides (MnO, Cr2O3, NiO, ...) into corresponding metals (Mn, Cr, Ni, ...) in place if extant.
If metals are expected as parts per million (ppm), set unitratio=10000 (default); if metals are as wt%, set unitratio = 1
See also `metalconversion`, c.f. `oxideconversion!`
"""
function metalconversion!(dataset::NamedTuple; unitratio::Number=10000)
# List of elements to convert
dest = (:Mn, :P, :Cr, :Ni, :Co, :Sr, :Ba, :Li, :S,)
source = (:MnO, :P2O5, :Cr2O3, :NiO, :CoO, :SrO, :BaO, :Li2O, :SO3)
conversionfactor = (0.7744619872751028, 0.4364268173666496, 0.6842080746199798, 0.785801615263874, 0.786486968277016, 0.8455993051534453, 0.8956541815613328, 0.46454031259412965, 0.4004646689233921)
# If source field exists, fill in destination from source
for i ∈ eachindex(source)
if haskey(dataset, source[i])
if haskey(dataset, dest[i]) # If destination field doesn't exist, make it.
metal, oxide = dataset[dest[i]], dataset[source[i]]
fillifnan!(metal, oxide, conversionfactor[i]*unitratio)
end
end
end
return dataset
end
metalconversion!(ds::Dict; kwargs...) = (metalconversion!(TupleDataset(ds); kwargs...); ds)
export metalconversion!
"""
```julia
carbonateconversion!(dataset::NamedTuple)
```
Convert carbonates (CaCO3, MgCO3) into corresponding metal oxides and CO2 if extant, in place,
as well as synchonizing TIC, TOC, TC, C and CO2. All are assumed to be reported in the same units,
(likely wt. %) except for C, which is assumed to be equivalent to unitratio * TC,
"""
function carbonateconversion!(ds::NamedTuple; unitratio=10000)
# Calculate CO2 if both CaCO3 and MgCO3 are reported
if haskey(ds, :CaCO3) && haskey(ds, :MgCO3) && haskey(ds, :CO2)
fillifnan!(ds.CO2, ds.CaCO3*0.43971009048182363 .+ ds.MgCO3*0.5219717006867268)
end
# Populate oxides and CO2 from carbonates and TIC
source = (:CaCO3, :CaCO3, :MgCO3, :MgCO3, :TIC,)
dest = (:CaO, :CO2, :MgO, :CO2, :CO2)
conversionfactor = (0.5602899095181764, 0.43971009048182363, 0.4780282993132732, 0.5219717006867268, 3.664057946882025)
for i in eachindex(source)
if haskey(ds, source[i])
if haskey(ds, dest[i])
d, s = ds[dest[i]], ds[source[i]]
fillifnan!(d, s, conversionfactor[i])
end
end
end
# Fill TC from C and TIC from CO2
if haskey(ds,:TC) && haskey(ds, :C)
fillifnan!(ds.TC, ds.C, 1e-4)
end
if haskey(ds,:TIC) && haskey(ds, :CO2)
fillifnan!(ds.TIC, ds.CO2, 0.27292144788565975)
end
# Synchronise TOC, TIC, TC
if haskey(ds, :TC) && haskey(ds, :TOC) && haskey(ds, :TIC)
fillifnan!(ds.TC, ds.TOC + ds.TIC)
fillifnan!(ds.TOC, ds.TC - ds.TIC)
nannegative!(ds.TOC)
fillifnan!(ds.TIC, ds.TC - ds.TOC)
nannegative!(ds.TIC)
if haskey(ds, :CO2)
# If we have new TIC values, fill CO2 again
fillifnan!(ds.CO2, ds.TIC, 3.664057946882025)
end
end
return ds
end
carbonateconversion!(ds::Dict) = (carbonateconversion!(TupleDataset(ds)); ds)
export carbonateconversion!
## --- Chemical Index of Alteration
# Chemical Index of Alteration as defined by Nesbitt and Young, 1982
# Note that CaO should be only igneous CaO excluding any Ca from calcite or apatite
function CIA(Al2O3::Number, CaO::Number, Na2O::Number, K2O::Number)
A = Al2O3 / 101.96007714
C = CaO / 56.0774
N = Na2O / 61.978538564
K = K2O / 94.19562
return A / (A + C + N + K) * 100
end
export CIA
# "Weathering Index of Parker" as defined by Parker, 1970
function WIP(Na2O::Number, MgO::Number, K2O::Number, CaO::Number)
Na = Na2O / 30.9895
Mg = MgO / 40.3044
K = K2O / 47.0980
Ca = CaO / 56.0774
# Denominator for each element is a measure of Nicholls' bond strengths
return (Na/0.35 + Mg/0.9 + K/0.25 + Ca/0.7) * 100
end
export WIP
## --- MELTS interface
"""
```julia
melts_configure(meltspath::String, scratchdir::String, composition::AbstractArray{<:Number},
\telements::AbstractArray{String},
\tT_range=(1400, 600),
\tP_range=(10000,10000);)
```
Configure and run a MELTS simulation using alphaMELTS.
Optional keyword arguments and defaults include:
batchstring::String = "1\nsc.melts\n10\n1\n3\n1\nliquid\n1\n1.0\n0\n10\n0\n4\n0\n"
A string defining the sequence of options that would be entered to produce
the desired calculation if running alphaMELTS at the command line. The
default string specifies a batch calculation starting at the liquidus.
dT = -10
The temperature step, in degrees, between each step of the MELTS calculation
dP = 0
The pressure step, in bar, between each step of the MELTS calculation
index = 1
An optional variable used to specify a unique suffix for the run directory name
version::String = "pMELTS"
A string specifying the desired version of MELTS. Options include `MELTS` and `pMELTS`.
mode::String = "isobaric"
A string specifying the desired calculation mode for MELTS. Options include
`isothermal`, `isobaric`, `isentropic`, `isenthalpic`, `isochoric`,
`geothermal` and `PTPath`.
fo2path::String = "FMQ"
A string specifying the oxygen fugacity buffer to follow, e.g., `FMQ` or `NNO+1`.
Available buffers include `IW`,`COH`,`FMQ`,`NNO`,`HM`, and `None`
fractionatesolids::Bool = false
Fractionate all solids? default is `false`
suppress::AbstractArray{String} = String[]
Supress individual phases (specify as strings in array, i.e. `["leucite"]`)
verbose::Bool = true
Print verbose MELTS output to terminal (else, write it to `melts.log`)
"""
function melts_configure(meltspath::String, scratchdir::String, composition::Collection{Number},
elements::AbstractArray{String}, T_range::Collection{Number}=(1400, 600), P_range::Collection{Number}=(10000,10000);
batchstring::String="1\nsc.melts\n10\n1\n3\n1\nliquid\n1\n1.0\n0\n10\n0\n4\n0\n",
dT=-10, dP=0, index=1, version="pMELTS",mode="isobaric",fo2path="FMQ",
fractionatesolids::Bool=false, suppress::AbstractArray{String}=String[], verbose::Bool=true)
############################ Default Settings ###############################
##MELTS or pMELTS
#version = "pMELTS"
##Mode (isothermal, isobaric, isentropic, isenthalpic, isochoric, geothermal or PTPath)
#mode = "isobaric"
## Set fO2 constraint, i.e. "IW","COH","FMQ","NNO","HM","None" as a string
#fo2path = "FMQ"
## Fractionate all solids? ("!" for no, "" for yes)
#fractionatesolids = "!"
# Mass retained during fractionation
massin = 0.001
# Ouptut temperatures in celcius? ("!" for no, "" for yes)
celciusoutput = ""
# Save all output? ("!" for no, "" for yes)
saveall = "!"
# Fractionate all water? ("!" for no, "" for yes)
fractionatewater = "!"
# Fractionate individual phases (specify as strings in cell array, i.e. {"olivine","spinel"})
fractionate = String[]
# Coninuous (fractional) melting? ("!" for no, "" for yes)
continuous = "!"
# Threshold above which melt is extracted (if fractionation is turned on)
minf = 0.005
# Do trace element calculations
dotrace = "!"
# Treat water as a trace element
dotraceh2o = "!"
# Initial trace compositionT
tsc = Float64[]
# Initial trace elements
telements = String[]
# Default global constraints
Pmax = 90000
Pmin = 2
Tmax = 3000
Tmin = 450
# Simulation number (for folder, etc)
########################## end Default Settings ############################
# Guess if intention is for calculation to end at Tf or Pf as a min or max
if last(T_range)<first(T_range)
Tmin=last(T_range)
end
if last(T_range)>first(T_range)
Tmax=last(T_range)
end
if last(P_range)<first(P_range)
Pmin=last(P_range)
end
if last(P_range)>first(P_range)
Pmax=last(P_range)
end
if fractionatesolids
fractionatesolids = ""
else
fractionatesolids = "!"
end
# Normalize starting composition
composition = composition./sum(composition)*100
# output prefixectory name
prefix = joinpath(scratchdir, "out$(index)/")
# Ensure directory exists and is empty
system("rm -rf $prefix; mkdir -p $prefix")
# Make .melts file containing the starting composition you want to run simulations on
fp = open(prefix*"sc.melts", "w")
for i ∈ eachindex(elements)
write(fp,"Initial Composition: $(elements[i]) $(trunc(composition[i],digits=4))\n")
end
for i ∈ eachindex(telements)
write(fp, "Initial Trace: $(telements[i]) $(trunc(tsc[i],digits=4))\n")
end
write(fp, "Initial Temperature: $(trunc(first(T_range),digits=2))\nInitial Pressure: $(trunc(first(P_range),digits=2))\nlog fo2 Path: $fo2path\n")
for i ∈ eachindex(fractionate)
write(fp,"Fractionate: $(fractionate[i])\n")
end
for i ∈ eachindex(suppress)
write(fp,"Suppress: $(suppress[i])\n")
end
close(fp)
# Make melts_env file to specify type of MELTS calculation
fp = open(prefix*"/melts_env.txt", "w")
write(fp, "! *************************************\n! Julia-generated environment file\n! *************************************\n\n" *
"! this variable chooses MELTS or pMELTS; for low-pressure use MELTS\n" *
"ALPHAMELTS_VERSION $version\n\n" *
"! do not use this unless fO2 anomalies at the solidus are a problem\n" *
"!ALPHAMELTS_ALTERNATIVE_FO2 true\n\n" *
"! use this if you want to buffer fO2 for isentropic, isenthalpic or isochoric mode\n! e.g. if you are doing isenthalpic AFC\n" *
"!ALPHAMELTS_IMPOSE_FO2 true\n\n" *
"! use if you want assimilation and fractional crystallization (AFC)\n" *
"!ALPHAMELTS_ASSIMILATE true\n\n" *
"! isothermal, isobaric, isentropic, isenthalpic, isochoric, geothermal or PTPath\n" *
"ALPHAMELTS_MODE $mode\n" *
"!ALPHAMELTS_PTPATH_FILE ptpath.txt\n\n" *
"! need to set DELTAP for polybaric paths; DELTAT for isobaric paths\nALPHAMELTS_DELTAP $(trunc(dP,digits=1))\n" *
"ALPHAMELTS_DELTAT $(trunc(dT,digits=1))\n" *
"ALPHAMELTS_MAXP $(trunc(Pmax,digits=1))\n" *
"ALPHAMELTS_MINP $(trunc(Pmin,digits=1))\n" *
"ALPHAMELTS_MAXT $(trunc(Tmax,digits=1))\n" *
"ALPHAMELTS_MINT $(trunc(Tmin,digits=1))\n\n" *
"! this one turns on fractional crystallization for all solids\n! use Fractionate: in the melts file instead for selective fractionation\n" *
"$(fractionatesolids)ALPHAMELTS_FRACTIONATE_SOLIDS true\n" *
"$(fractionatesolids)ALPHAMELTS_MASSIN $massin\n\n" *
"! free water is unlikely but can be extracted\n" *
"$(fractionatewater)ALPHAMELTS_FRACTIONATE_WATER true\n" *
"$(fractionatewater)ALPHAMELTS_MINW 0.005\n\n" *
"! the next one gives an output file that is always updated, even for single calculations\n" *
"$(saveall)ALPHAMELTS_SAVE_ALL true\n" *
"!ALPHAMELTS_SKIP_FAILURE true\n\n" *
"! this option converts the output temperature to celcius, like the input\n" *
"$(celciusoutput)ALPHAMELTS_CELSIUS_OUTPUT true\n\n" *
"! the next two turn on and off fractional melting\n" *
"$(continuous)ALPHAMELTS_CONTINUOUS_MELTING true\n" *
"$(continuous)ALPHAMELTS_MINF $minf\n" *
"$(continuous)ALPHAMELTS_INTEGRATE_FILE integrate.txt\n\n" *
"! the next two options refer to the trace element engine\n" *
"$(dotrace)ALPHAMELTS_DO_TRACE true\n" *
"$(dotraceh2o)ALPHAMELTS_DO_TRACE_H2O true\n")
close(fp)
# Make a batch file to run the above .melts file starting from the liquidus
fp = open(prefix*"/batch.txt", "w")
write(fp,batchstring)
close(fp)
# Run the command
# Edit the following line(s to make sure you have a correct path to the "run_alphamelts.command" perl script
if verbose
system("cd " * prefix * "; " * meltspath * " -f melts_env.txt -b batch.txt")
else
system("cd " * prefix * "; " * meltspath * " -f melts_env.txt -b batch.txt &>./melts.log")
end
return 0
end
export melts_configure
"""
```julia
melts_query(scratchdir::String; index=1)
```
Read all phase proportions from `Phase_main_tbl.txt` in specified MELTS run directory
Returns an elementified dictionary
"""
function melts_query(scratchdir::String; index=1, importas=:Dict)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
if importas==:Dict
melts = Dict{String, Union{Vector{String}, Dict}}()
else
melts = Dict{String, Union{Vector{String}, NamedTuple}}()
end
if isfile(prefix*"/Phase_main_tbl.txt")
data = readdlm(prefix*"/Phase_main_tbl.txt", ' ', skipblanks=false)
pos = findall(all(isempty.(data), dims=2) |> vec)
melts["minerals"] = Array{String}(undef, length(pos)-1)
for i=1:(length(pos)-1)
name = data[pos[i]+1,1]
melts[name] = elementify(data[pos[i]+2:pos[i+1]-1,:], skipnameless=true, importas=importas)
melts["minerals"][i] = name
end
end
return melts
end
export melts_query
"""
```julia
melts_query_modes(scratchdir::String; index=1)
```
Read modal phase proportions from `Phase_mass_tbl.txt` in specified MELTS run
Returns an elementified dictionary
"""
function melts_query_modes(scratchdir::String; index=1, importas=:Dict)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
# Read results and return them if possible
if isfile(prefix*"/Phase_mass_tbl.txt")
# Read data as an Array{Any}
data = readdlm(prefix*"Phase_mass_tbl.txt", ' ', skipstart=1)
# Convert to a dictionary
data = elementify(data, standardize=true, skipnameless=true, importas=importas)
else
# Return empty dictionary if file doesn't exist
data = importas==:Dict ? Dict() : ()
end
return data
end
export melts_query_modes
"""
```julia
melts_clean_modes(scratchdir::String; index=1)
```
Read and parse / clean-up modal phase proportions from specified MELTS run directory
Returns an elementified dictionary
"""
function melts_clean_modes(scratchdir::String; index=1)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
# Read results and return them if possible
if isfile(prefix*"/Phase_mass_tbl.txt")
# Read data as an Array{Any}
data = readdlm(prefix*"Phase_mass_tbl.txt", ' ', skipstart=1)
# Convert to a dictionary
data = elementify(data, standardize=true, skipnameless=true, importas=:Dict)
# Start by transferring over all the non-redundant elements
modes = typeof(data)()
for e in data["elements"]
m = replace(e, r"_.*" => s"")
if haskey(modes, m)
modes[m] .+= data[e]
else
modes[m] = copy(data[e])
end
end
# Add the sum of all solids
modes["solids"] = zeros(size(data["Temperature"]))
for e in data["elements"][4:end]
if !contains(e, "water") && !contains(e, "liquid")
modes["solids"] .+= data[e]
end
end
# Get full mineral compositions, add feldspar and oxides
melts = melts_query(scratchdir, index=index)
if containsi(melts["minerals"],"feldspar")
modes["anorthite"] = zeros(size(modes["Temperature"]))
modes["albite"] = zeros(size(modes["Temperature"]))
modes["orthoclase"] = zeros(size(modes["Temperature"]))
end
An_Ca = (238.12507+40.0784) / (15.999+40.0784)
Ab_Na = (239.22853+22.98977*2) / (15.999+22.98977*2)
Or_K = (239.22853+39.09831*2) / (15.999+39.09831*2)
if containsi(melts["minerals"],"rhm_oxide")
modes["ilmenite"] = zeros(size(modes["Temperature"]))
modes["magnetite"] = zeros(size(modes["Temperature"]))
modes["hematite"] = zeros(size(modes["Temperature"]))
end
for m in melts["minerals"]
if containsi(m,"feldspar")
t = vec(findclosest(melts[m]["Temperature"],modes["Temperature"]))
AnAbOr = [melts[m]["CaO"]*An_Ca melts[m]["Na2O"]*Ab_Na melts[m]["K2O"]*Or_K] |> x -> x ./ sum(x, dims=2)
modes["anorthite"][t] .+= AnAbOr[:,1] .* melts[m]["mass"]
modes["albite"][t] .+= AnAbOr[:,2] .* melts[m]["mass"]
modes["orthoclase"][t] .+= AnAbOr[:,3] .* melts[m]["mass"]
elseif containsi(m,"rhm_oxide")
t = vec(findclosest(melts[m]["Temperature"],modes["Temperature"]))
Ilmenite = Vector{Float64}(undef, length(t))
Magnetite = Vector{Float64}(undef, length(t))
if haskey(melts[m],"MnO")
Ilmenite .= (melts[m]["TiO2"] + melts[m]["MnO"]+(melts[m]["TiO2"]*(71.8444/79.8768) - melts[m]["MnO"]*(71.8444/70.9374))) / 100
Magnetite .= (melts[m]["FeO"] - (melts[m]["TiO2"])*71.8444/79.8768) * (1+159.6882/71.8444)/100
else
Ilmenite .= (melts[m]["TiO2"] + melts[m]["TiO2"]*71.8444/79.8768) / 100
Magnetite .= (melts[m]["FeO"] - melts[m]["TiO2"]*71.8444/79.8768) * (1+159.6882/71.8444)/100
end
Magnetite[Magnetite.<0] .= 0
Hematite = (melts[m]["Fe2O3"] - Magnetite*100*159.6882/231.5326)/100
modes["ilmenite"][t] .+= melts[m]["mass"] .* Ilmenite
modes["magnetite"][t] .+= melts[m]["mass"] .* Magnetite
modes["hematite"][t] .+= melts[m]["mass"] .* Hematite
end
end
minerals = sort(collect(keys(modes)))
modes["elements"] = ["Pressure","Temperature","mass","solids","liquid"] ∪ minerals[.!containsi.(minerals, "feldspar") .& .!containsi.(minerals, "rhm")]
else
# Return empty dictionary if file doesn't exist
modes = Dict()
end
return modes
end
export melts_clean_modes
"""
```julia
melts_query_liquid(scratchdir::String; index=1)
```
Read liquid composition from `Liquid_comp_tbl.txt` in specified MELTS run directory
Returns an elementified dictionary
"""
function melts_query_liquid(scratchdir::String; index=1, importas=:Dict)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
# Read results and return them if possible
if isfile(prefix*"/Liquid_comp_tbl.txt")
# Read data as an Array{Any}
data = readdlm(prefix*"Liquid_comp_tbl.txt", ' ', skipstart=1)
# Convert to a dictionary
data = elementify(data, standardize=true, skipnameless=true, importas=importas)
else
# Return empty dictionary if file doesn't exist
data = importas==:Dict ? Dict() : ()
end
return data
end
export melts_query_liquid
"""
```julia
melts_query_solid(scratchdir::String; index=1)
```
Read solid composition from `Solid_comp_tbl.txt` in specified MELTS run directory
Returns an elementified dictionary
"""
function melts_query_solid(scratchdir::String; index=1, importas=:Dict)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
# Read results and return them if possible
if isfile(prefix*"/Solid_comp_tbl.txt")
# Read data as an Array{Any}
data = readdlm(prefix*"Solid_comp_tbl.txt", ' ', skipstart=1)
# Convert to a dictionary
data = elementify(data, standardize=true, skipnameless=true, importas=importas)
else
# Return empty dictionary if file doesn't exist
data = importas==:Dict ? Dict() : ()
end
return data
end
export melts_query_solid
"""
```julia
melts_query_system(scratchdir::String; index=1, importas=:Dict)
```
Read system thermodynamic and composition data from `System_main_tbl.txt` in
specified MELTS run directory. Returns an elementified dictionary or tuple.
"""
function melts_query_system(scratchdir::String; index=1, importas=:Dict)
prefix = joinpath(scratchdir, "out$(index)/") # path to data files
# Read results and return them if possible
if isfile(prefix*"/System_main_tbl.txt")
# Read data as an Array{Any}
data = readdlm(prefix*"System_main_tbl.txt", ' ', skipstart=1)
# Convert to a dictionary
data = elementify(data, standardize=true, skipnameless=true, importas=importas)
else
# Return empty dictionary if file doesn't exist
data = importas==:Dict ? Dict() : ()
end
return data
end
export melts_query_system
## -- Perplex interface: 1. Configuration
"""
```julia
perplex_configure_geotherm(perplexdir::String, scratchdir::String, composition::Collection{Number},
\telements::String=["SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"],
\tP_range=(280,28000), T_surf::Number=273.15, geotherm::Number=0.1;
\tdataset::String="hp02ver.dat",
\tindex::Integer=1,
\tnpoints::Integer=100,
\tsolution_phases::String="O(HP)\\nOpx(HP)\\nOmph(GHP)\\nGt(HP)\\noAmph(DP)\\ncAmph(DP)\\nT\\nB\\nChl(HP)\\nBio(TCC)\\nMica(CF)\\nCtd(HP)\\nIlHm(A)\\nSp(HP)\\nSapp(HP)\\nSt(HP)\\nfeldspar_B\\nDo(HP)\\nF\\n",
\texcludes::String="ts\\nparg\\ngl\\nged\\nfanth\\ng\\n",
\tmode_basis::String="vol", #["vol", "wt", "mol"]
\tcomposition_basis::String="wt", #["vol", "wt", "mol"]
\tfluid_eos::Integer=5)
```
Set up a PerpleX calculation for a single bulk composition along a specified
geothermal gradient and pressure (depth) range. P specified in bar and T_surf
in Kelvin, with geothermal gradient in units of Kelvin/bar
"""
function perplex_configure_geotherm(perplexdir::String, scratchdir::String, composition::Collection{Number},
elements::Collection{String}=["SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"],
P_range::NTuple{2,Number}=(280,28000), T_surf::Number=273.15, geotherm::Number=0.1;
dataset::String="hp02ver.dat",
index::Integer=1,
npoints::Integer=100,
solution_phases::String="O(HP)\nOpx(HP)\nOmph(GHP)\nGt(HP)\noAmph(DP)\ncAmph(DP)\nT\nB\nChl(HP)\nBio(TCC)\nMica(CF)\nCtd(HP)\nIlHm(A)\nSp(HP)\nSapp(HP)\nSt(HP)\nfeldspar_B\nDo(HP)\nF\n",
excludes::String="ts\nparg\ngl\nged\nfanth\ng\n",
mode_basis::String="vol",
composition_basis::String="wt",
fluid_eos::Integer=5
)
build = joinpath(perplexdir, "build")# path to PerpleX build
vertex = joinpath(perplexdir, "vertex")# path to PerpleX vertex
#Configure working directory
prefix = joinpath(scratchdir, "out$(index)/")
system("rm -rf $prefix; mkdir -p $prefix")
# Place required data files
system("cp $(joinpath(perplexdir,dataset)) $prefix")
system("cp $(joinpath(perplexdir,"perplex_option.dat")) $prefix")
system("cp $(joinpath(perplexdir,"solution_model.dat")) $prefix")
# Edit perplex_option.dat to specify number of nodes at which to solve
system("sed -e \"s/1d_path .*|/1d_path $npoints $npoints |/\" -i.backup $(prefix)perplex_option.dat")
# Edit perplex_option.dat to output all seismic properties
#println("editing perplex options ")
system("sed -e \"s/seismic_output .*|/seismic_output all |/\" -i.backup $(prefix)perplex_option.dat")
# Specify whether we want volume or weight percentages
system("sed -e \"s/proportions .*|/proportions $mode_basis |/\" -i.backup $(prefix)perplex_option.dat")
system("sed -e \"s/composition_system .*|/composition_system $composition_basis |/\" -i.backup $(prefix)perplex_option.dat")
system("sed -e \"s/composition_phase .*|/composition_phase $composition_basis |/\" -i.backup $(prefix)perplex_option.dat")
# Create build batch file.
fp = open(prefix*"build.bat", "w")
# Name, components, and basic options. P-T conditions.
# default fluid_eos = 5: Holland and Powell (1998) "CORK" fluid equation of state
elementstring = join(elements .* "\n")
write(fp,"$index\n$dataset\nperplex_option.dat\nn\n3\nn\nn\nn\n$elementstring\n$fluid_eos\nn\ny\n2\n1\n$T_surf\n$geotherm\n$(first(P_range))\n$(last(P_range))\ny\n") # v6.8.7
# write(fp,"$index\n$dataset\nperplex_option.dat\nn\nn\nn\nn\n$elementstring\n5\n3\nn\ny\n2\n1\n$T_surf\n$geotherm\n$(first(P_range))\n$(last(P_range))\ny\n") # v6.8.1
# Whole-rock composition
for i ∈ eachindex(composition)
write(fp,"$(composition[i]) ")
end
# Solution model
if length(excludes) > 0
write(fp,"\nn\ny\nn\n$excludes\ny\nsolution_model.dat\n$solution_phases\nGeothermal")
else
write(fp,"\nn\nn\ny\nsolution_model.dat\n$(solution_phases)\nGeothermal")
end
close(fp)
# build PerpleX problem definition
system("cd $prefix; $build < build.bat > build.log")
println("Built problem definition")
# Run PerpleX vertex calculations
result = system("cd $prefix; echo $index | $vertex > vertex.log")
return result
end
export perplex_configure_geotherm
"""
```julia
perplex_configure_isobar(perplexdir::String, scratchdir::String, composition::Collection{Number},
\telements::String=["SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"]
\tP::Number=10000, T_range::NTuple{2,Number}=(500+273.15, 1500+273.15);
\tdataset::String="hp11ver.dat",
\tindex::Integer=1,
\tnpoints::Integer=100,
\tsolution_phases::String="O(HP)\\nOpx(HP)\\nOmph(GHP)\\nGt(HP)\\noAmph(DP)\\ncAmph(DP)\\nT\\nB\\nChl(HP)\\nBio(TCC)\\nMica(CF)\\nCtd(HP)\\nIlHm(A)\\nSp(HP)\\nSapp(HP)\\nSt(HP)\\nfeldspar_B\\nDo(HP)\\nF\\n",
\texcludes::String="ts\\nparg\\ngl\\nged\\nfanth\\ng\\n",
\tmode_basis::String="wt", #["vol", "wt", "mol"]
\tcomposition_basis::String="wt", #["vol", "wt", "mol"]
\tnonlinear_subdivision::Bool=false,
\tfluid_eos::Integer=5)
```
Set up a PerpleX calculation for a single bulk composition along a specified
isobaric temperature gradient. P specified in bar and T_range in Kelvin
"""
function perplex_configure_isobar(perplexdir::String, scratchdir::String, composition::Collection{Number},
elements::Collection{String}=("SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"),
P::Number=10000, T_range::NTuple{2,Number}=(500+273.15, 1500+273.15);
dataset::String="hp11ver.dat",
index::Integer=1,
npoints::Integer=100,
solution_phases::String="O(HP)\nOpx(HP)\nOmph(GHP)\nGt(HP)\noAmph(DP)\ncAmph(DP)\nT\nB\nChl(HP)\nBio(TCC)\nMica(CF)\nCtd(HP)\nIlHm(A)\nSp(HP)\nSapp(HP)\nSt(HP)\nfeldspar_B\nDo(HP)\nF\n",
excludes::String="ts\nparg\ngl\nged\nfanth\ng\n",
mode_basis::String="wt",
composition_basis::String="wt",
nonlinear_subdivision::Bool=false,
fluid_eos::Integer=5
)
build = joinpath(perplexdir, "build")# path to PerpleX build
vertex = joinpath(perplexdir, "vertex")# path to PerpleX vertex
#Configure working directory
prefix = joinpath(scratchdir, "out$(index)/")
system("rm -rf $prefix; mkdir -p $prefix")
# Place required data files
system("cp $(joinpath(perplexdir,dataset)) $prefix")
system("cp $(joinpath(perplexdir,"perplex_option.dat")) $prefix")
system("cp $(joinpath(perplexdir,"solution_model.dat")) $prefix")
# Edit perplex_option.dat to specify number of nodes at which to solve
system("sed -e \"s/1d_path .*|/1d_path $npoints $npoints |/\" -i.backup $(prefix)perplex_option.dat")
# Specify whether we want volume or weight percentages
system("sed -e \"s/proportions .*|/proportions $mode_basis |/\" -i.backup $(prefix)perplex_option.dat")
system("sed -e \"s/composition_system .*|/composition_system $composition_basis |/\" -i.backup $(prefix)perplex_option.dat")
system("sed -e \"s/composition_phase .*|/composition_phase $composition_basis |/\" -i.backup $(prefix)perplex_option.dat")
# Turn on nonlinear subdivision and change resolution
if nonlinear_subdivision
system("sed -e \"s/non_linear_switch .*|/non_linear_switch T |/\" -i.backup $(prefix)perplex_option.dat")
system("sed -e \"s:initial_resolution .*|:initial_resolution 1/2 1/4 |:\" -i.backup $(prefix)perplex_option.dat")
end
# Create build batch file
# Options based on Perplex v6.8.7
fp = open(prefix*"build.bat", "w")
# Name, components, and basic options. P-T conditions.
# default fluid_eos = 5: Holland and Powell (1998) "CORK" fluid equation of state
elementstring = join(elements .* "\n")
write(fp,"$index\n$dataset\nperplex_option.dat\nn\n3\nn\nn\nn\n$elementstring\n$fluid_eos\nn\nn\n2\n$(first(T_range))\n$(last(T_range))\n$P\ny\n") # v6.8.7
# write(fp,"$index\n$dataset\nperplex_option.dat\nn\nn\nn\nn\n$elementstring\n$fluid_eos\n3\nn\nn\n2\n$(first(T_range))\n$(last(T_range))\n$P\ny\n") # v6.8.1
# Whole-rock composition
for i ∈ eachindex(composition)
write(fp,"$(composition[i]) ")
end
# Solution model
write(fp,"\nn\ny\nn\n$excludes\ny\nsolution_model.dat\n$solution_phases\nIsobaric")
close(fp)
# build PerpleX problem definition
system("cd $prefix; $build < build.bat > build.log")
# Run PerpleX vertex calculations
result = system("cd $prefix; printf \"$index\ny\ny\ny\ny\ny\ny\ny\ny\n\" | $vertex > vertex.log")
return result
end
export perplex_configure_isobar
"""
```julia
perplex_configure_path(perplexdir::String, scratchdir::String, composition::Collection{Number}, PTdir::String="",
\telements::String=("SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"),
\tT_range::NTuple{2,Number}=(500+273.15, 1050+273.15);
\tdataset::String="hp11ver.dat",
\tindex::Integer=1,
\tsolution_phases::String="O(HP)\\nOpx(HP)\\nOmph(GHP)\\nGt(HP)\\noAmph(DP)\\ncAmph(DP)\\nT\\nB\\nChl(HP)\\nBio(TCC)\\nMica(CF)\\nCtd(HP)\\nIlHm(A)\\nSp(HP)\\nSapp(HP)\\nSt(HP)\\nfeldspar_B\\nDo(HP)\\nF\\n",
\texcludes::String="ts\\nparg\\ngl\\nged\\nfanth\\ng\\n",
\tmode_basis::String="wt", #["vol", "wt", "mol"]
\tcomposition_basis::String="wt", #["vol", "wt", "mol"]
\tnonlinear_subdivision::Bool=false,
\tfluid_eos::Integer=5)
```
Set up a PerpleX calculation for a single bulk composition along a specified
pressure–temperature path with T as the independent variable.
P specified in bar and T_range in Kelvin
"""
function perplex_configure_path(perplexdir::String, scratchdir::String, composition::Collection{Number}, PTdir::String="",
elements::Collection{String}=("SIO2","TIO2","AL2O3","FEO","MGO","CAO","NA2O","K2O","H2O"),
T_range::NTuple{2,Number}=(500+273.15, 1050+273.15);