-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathatmosphere.c
1204 lines (918 loc) · 43.7 KB
/
atmosphere.c
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
#include <petsc.h>
#include "atmosphere.h"
#include "ctx.h"
#include "dimensionalisablefield.h"
#include "reaction.h"
#include "util.h"
static PetscScalar get_grey_body_flux( const Atmosphere *, const AtmosphereParameters, const FundamentalConstants );
static PetscScalar get_steam_atmosphere_zahnle_1988_flux( const Atmosphere *, const ScalingConstants );
static PetscScalar get_emissivity_abe_matsui( Atmosphere *, const AtmosphereParameters);
static PetscScalar get_emissivity_from_flux( const Atmosphere *, const AtmosphereParameters, const FundamentalConstants, PetscScalar );
static PetscErrorCode set_total_surface_pressure( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_total_surface_pressure_time_derivative( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_volume_mixing_ratios( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_volatile_masses_in_atmosphere( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_volatile_masses_in_liquid( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_volatile_masses_in_solid( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_volatile_masses_reactions( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_jeans_escape( Atmosphere *, const AtmosphereParameters, const FundamentalConstants );
static PetscErrorCode set_zahnle_escape( Atmosphere *, const AtmosphereParameters );
static PetscScalar get_pressure_dependent_kabs( const AtmosphereParameters, PetscInt );
static PetscErrorCode set_jeans( Atmosphere *, const AtmosphereParameters, const FundamentalConstants, PetscInt );
static PetscErrorCode set_column_density_volatile( Atmosphere *, const AtmosphereParameters, const FundamentalConstants, PetscInt );
static PetscErrorCode set_Knudsen_number( Atmosphere *, const AtmosphereParameters, PetscInt );
static PetscErrorCode set_R_thermal_escape( Atmosphere *, const AtmosphereParameters, PetscInt );
static PetscErrorCode set_f_thermal_escape( Atmosphere *, PetscInt );
static PetscScalar get_x_from_solubility_power_law( PetscScalar partialp, PetscScalar henry, PetscScalar henry_pow );
static PetscErrorCode JSON_add_volatile( DM, Parameters const, VolatileParameters const, Volatile const *, Atmosphere const *, char const *name, cJSON * );
static PetscErrorCode JSON_add_atm_struct( Atmosphere *, const AtmosphereParameters, const FundamentalConstants, cJSON * );
static PetscErrorCode set_atm_struct_tau( Atmosphere * );
static PetscErrorCode set_atm_struct_temp( Atmosphere *, const AtmosphereParameters, const FundamentalConstants );
static PetscErrorCode set_atm_struct_pressure( Atmosphere *, const AtmosphereParameters );
static PetscErrorCode set_atm_struct_depth( Atmosphere *, const AtmosphereParameters, const FundamentalConstants );
PetscErrorCode initialise_atmosphere( Atmosphere *A, const AtmosphereParameters Ap, const ScalingConstants SC )
{
PetscErrorCode ierr;
PetscScalar scaling;
PetscFunctionBeginUser;
// pointer to da, so we can easily access it within the atmosphere structure
const PetscInt stencilWidth = 1;
const PetscInt dof = 1;
const PetscInt numpts = 500; // hard-coded for atmosphere structure output
ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,numpts,dof,stencilWidth,NULL,&A->da_atm);CHKERRQ(ierr);
ierr = DMSetUp(A->da_atm);CHKERRQ(ierr);
/* create dimensionalisable fields for outputting atmosphere structure */
/* this does not really need to be a DimensionalisableField, and PS
proposed creating a new structure call a DimenisonalisableValue instead */
scaling = 1.0;
ierr = DimensionalisableFieldCreate(&A->atm_struct[0],A->da_atm,&scaling,PETSC_FALSE);CHKERRQ(ierr);
ierr = DimensionalisableFieldGetGlobalVec(A->atm_struct[0],&A->atm_struct_tau);
ierr = DimensionalisableFieldSetName(A->atm_struct[0],"atm_struct_tau");CHKERRQ(ierr);
ierr = DimensionalisableFieldSetUnits(A->atm_struct[0],"None");CHKERRQ(ierr);
scaling = SC->TEMP;
ierr = DimensionalisableFieldCreate(&A->atm_struct[1],A->da_atm,&scaling,PETSC_FALSE);CHKERRQ(ierr);
ierr = DimensionalisableFieldGetGlobalVec(A->atm_struct[1],&A->atm_struct_temp);
ierr = DimensionalisableFieldSetName(A->atm_struct[1],"atm_struct_temp");CHKERRQ(ierr);
ierr = DimensionalisableFieldSetUnits(A->atm_struct[1],"K");CHKERRQ(ierr);
scaling = SC->PRESSURE / 1.0E5; // bar
ierr = DimensionalisableFieldCreate(&A->atm_struct[2],A->da_atm,&scaling,PETSC_FALSE);CHKERRQ(ierr);
ierr = DimensionalisableFieldGetGlobalVec(A->atm_struct[2],&A->atm_struct_pressure);
ierr = DimensionalisableFieldSetName(A->atm_struct[2],"atm_struct_pressure");CHKERRQ(ierr);
ierr = DimensionalisableFieldSetUnits(A->atm_struct[2],"bar");CHKERRQ(ierr);
scaling = SC->RADIUS;
ierr = DimensionalisableFieldCreate(&A->atm_struct[3],A->da_atm,&scaling,PETSC_FALSE);CHKERRQ(ierr);
ierr = DimensionalisableFieldGetGlobalVec(A->atm_struct[3],&A->atm_struct_depth);
ierr = DimensionalisableFieldSetName(A->atm_struct[3],"atm_struct_depth");CHKERRQ(ierr);
ierr = DimensionalisableFieldSetUnits(A->atm_struct[3],"m");CHKERRQ(ierr);
PetscFunctionReturn(0);
}
PetscErrorCode destroy_atmosphere( Atmosphere *A )
{
PetscErrorCode ierr;
PetscInt i;
PetscFunctionBeginUser;
for (i=0;i<NUMATMSTRUCTVECS;++i){
ierr = DimensionalisableFieldDestroy(&A->atm_struct[i]);CHKERRQ(ierr);
}
ierr = DMDestroy(&A->da_atm);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_atm_struct_tau( Atmosphere *A )
{
/* builds an evenly-spaced profile of optical depth from unity at
the top to the surface value */
PetscErrorCode ierr;
PetscScalar const logtau_min = -6; // hard-coded here
PetscScalar const logtau_max = PetscLog10Real(A->tau); // surface optical depth
PetscScalar tau,logdtau;
PetscInt i,ilo,w,ihi,numpts;
PetscFunctionBeginUser;
ierr = DMDAGetInfo(A->da_atm,NULL,&numpts,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
ierr = DMDAGetCorners(A->da_atm,&ilo,0,0,&w,0,0);CHKERRQ(ierr);
ihi = ilo + w;
logdtau = (logtau_max - logtau_min) / (numpts-1);
for(i=ilo; i<ihi; ++i){
tau = PetscPowScalar(10.0,logtau_min + i*logdtau);
ierr = VecSetValues( A->atm_struct_tau, 1, &i, &tau, INSERT_VALUES );CHKERRQ(ierr);
}
ierr = VecAssemblyBegin( A->atm_struct_tau );CHKERRQ(ierr);
ierr = VecAssemblyEnd( A->atm_struct_tau );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_atm_struct_temp( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC )
{
PetscErrorCode ierr;
PetscScalar TO4, Teq4;
PetscFunctionBeginUser;
TO4 = A->Fatm / FC->STEFAN_BOLTZMANN;
Teq4 = PetscPowScalar(Ap->teqm,4.0);
ierr = VecCopy( A->atm_struct_tau, A->atm_struct_temp ); CHKERRQ(ierr);
ierr = VecShift( A->atm_struct_temp, 1.0 ); CHKERRQ(ierr);
ierr = VecScale( A->atm_struct_temp, 0.5 ); CHKERRQ(ierr);
ierr = VecScale( A->atm_struct_temp, TO4 ); CHKERRQ(ierr);
ierr = VecShift( A->atm_struct_temp, Teq4 ); CHKERRQ(ierr);
ierr = VecPow( A->atm_struct_temp, 1.0/4.0 ); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_atm_struct_pressure( Atmosphere *A, const AtmosphereParameters Ap )
{
PetscErrorCode ierr;
PetscInt i;
PetscScalar kabs_tot,kabs;
PetscFunctionBeginUser;
/* effective absorption coefficient */
kabs_tot = 0.0;
for( i=0; i<Ap->n_volatiles; ++i) {
kabs = get_pressure_dependent_kabs( Ap, i );
kabs_tot += A->volatiles[i].mixing_ratio * kabs;
}
ierr = VecCopy( A->atm_struct_tau, A->atm_struct_pressure ); CHKERRQ(ierr);
ierr = VecScale( A->atm_struct_pressure, -(*Ap->gravity_ptr) ); CHKERRQ(ierr); // note negative gravity
ierr = VecScale( A->atm_struct_pressure, 2.0 ); CHKERRQ(ierr);
ierr = VecScale( A->atm_struct_pressure, 1.0/3.0 ); CHKERRQ(ierr);
ierr = VecScale( A->atm_struct_pressure, 1.0/kabs_tot); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_jeans( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, PetscInt i )
{
/* surface Jeans parameter */
Volatile *V = &A->volatiles[i];
VolatileParameters const Vp = Ap->volatile_parameters[i];
PetscFunctionBeginUser;
if(Vp->jeans_value < 0.0){
V->jeans = -(*Ap->gravity_ptr) * (*Ap->radius_ptr) * (Vp->molar_mass/FC->AVOGADRO); // note negative gravity
V->jeans /= FC->BOLTZMANN * A->tsurf;
}
else{
V->jeans = Vp->jeans_value;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_column_density_volatile( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, PetscInt i )
{
/* see Johnson et al. (2015), Astrophys. J. */
PetscFunctionBeginUser;
Volatile *V = &A->volatiles[i];
VolatileParameters const Vp = Ap->volatile_parameters[i];
V->column_density = V->p;
V->column_density /= -(*Ap->gravity_ptr) * (Vp->molar_mass/FC->AVOGADRO);
PetscFunctionReturn(0);
}
static PetscErrorCode set_Knudsen_number( Atmosphere *A, const AtmosphereParameters Ap, PetscInt i )
{
/* see Johnson et al. (2015), Astrophys. J. */
Volatile *V = &A->volatiles[i];
VolatileParameters const Vp = Ap->volatile_parameters[i];
PetscFunctionBeginUser;
V->Knudsen = V->jeans * Vp->cross_section * V->column_density;
V->Knudsen = 1.0 / V->Knudsen;
PetscFunctionReturn(0);
}
static PetscErrorCode set_R_thermal_escape( Atmosphere *A, const AtmosphereParameters Ap, PetscInt i )
{
/* see Johnson et al. (2015), Astrophys. J. */
PetscScalar R1, R2, Rfit;
Volatile *V = &A->volatiles[i];
VolatileParameters const Vp = Ap->volatile_parameters[i];
PetscFunctionBeginUser;
if(Vp->R_thermal_escape_value < 0.0 ){
R1 = PetscPowScalar( 1.0 / V->Knudsen, 0.09 );
R2 = 70 * (V->Knudsen * PetscExpReal(V->jeans)) / PetscPowScalar(V->jeans,2.55);
Rfit = 1.0 / R1 + 1.0 / R2;
Rfit = 1.0 / Rfit;
V->R_thermal_escape = Rfit;
}
else{
V->R_thermal_escape = Vp->R_thermal_escape_value;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_f_thermal_escape( Atmosphere *A, PetscInt i )
{
/* thermal escape prefactor for atmospheric growth rate */
/* see Johnson et al. (2015), Astrophys. J. */
Volatile *V = &A->volatiles[i];
PetscFunctionBeginUser;
V->f_thermal_escape = 1.0 + V->R_thermal_escape * (1.0+V->jeans) * PetscExpReal(-V->jeans);
/* As jeans-->infty the limit looks OK, but what about
as jeans-->0? In reality, this would denote a switch to hydrodynamic
escape */
PetscFunctionReturn(0);
}
static PetscErrorCode set_total_surface_pressure( Atmosphere *A, const AtmosphereParameters Ap )
{
PetscInt i;
Volatile *V;
PetscFunctionBeginUser;
/* total surface pressure */
A->psurf = 0.0;
for (i=0; i<Ap->n_volatiles; ++i) {
V = &A->volatiles[i];
/* partial pressure of volatile */
A->psurf += V->p;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_total_surface_pressure_time_derivative( Atmosphere *A, const AtmosphereParameters Ap )
{
/* dpsurf/dt */
PetscInt k;
PetscFunctionBeginUser;
/* total surface pressure time derivative */
A->dpsurfdt = 0.0;
for (k=0; k<Ap->n_volatiles; ++k) {
A->dpsurfdt += A->volatiles[k].dpdt;
}
PetscFunctionReturn(0);
}
static PetscScalar get_x_from_solubility_power_law( PetscScalar partialp, PetscScalar henry, PetscScalar henry_pow )
{
/* Solubility power law (Henry-like with exponent) */
return henry * PetscPowScalar( partialp, 1.0/henry_pow );
}
static PetscScalar get_dxdp_from_solubility_power_law( PetscScalar partialp, PetscScalar henry, PetscScalar henry_pow )
{
/* dxdp from solubility power law (Henry-like with exponent) */
return (henry / henry_pow ) * PetscPowScalar( partialp, 1.0/henry_pow - 1.0 );
}
PetscErrorCode set_volatile_abundances_from_partial_pressure( Atmosphere *A, const AtmosphereParameters Ap, const ScalingConstants SC )
{
/* This function contains the solubility laws. For each solubility law, you must give the relationship
between x and p. NOTE: for every solubility law, you must also specify dx/dp (hence dx/dt) in
get_dpdt() */
PetscInt i;
Volatile *V;
VolatileParameters Vp;
PetscScalar pbar=0, Tkel=0;
PetscFunctionBeginUser;
for (i=0; i<Ap->n_volatiles; ++i) {
V = &A->volatiles[i];
Vp = Ap->volatile_parameters[i];
switch( Vp->SOLUBILITY ){
case 1:
/* Modified Henry's law (default) */
/* x = henry * partialp ** (1/beta) */
V->x = get_x_from_solubility_power_law( V->p, Vp->henry, Vp->henry_pow );
break;
case 2:
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported SOLUBILITY value %d provided",Vp->SOLUBILITY);
break;
case 3:
/* CO2 for Basalt from Dixon et al. (1995) */
/* need T in K and P in bar for this solubility formulation */
pbar = V->p * SC->PRESSURE * 1.0E-5; // bar
Tkel = A->tsurf * SC->TEMP; // Kelvin
V->x = (3.8E-7)*pbar*PetscExpScalar( -23*(pbar-1)/(83.15*Tkel) );
V->x = 1.0E4 * (4400 * V->x ) / (36.6 - 44 * V->x ); // ppm
V->x *= 1.0E-6 / SC->VOLATILE;
break;
/* add more cases to include more solubility laws */
/* could also interface with self-consistent solubility calculation (PERPLEX) */
default:
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported SOLUBILITY value %d provided",Vp->SOLUBILITY);
}
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_volatile_masses_in_atmosphere( Atmosphere *A, const AtmosphereParameters Ap )
{
/* mass of volatiles in atmosphere */
PetscInt i;
PetscFunctionBeginUser;
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].mass_atmos = PetscSqr((*Ap->radius_ptr)) * A->volatiles[i].p / -(*Ap->gravity_ptr);
A->volatiles[i].mass_atmos /= (*Ap->VOLATILE_ptr);
A->volatiles[i].mass_atmos *= Ap->volatile_parameters[i]->molar_mass / A->molar_mass;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_volatile_masses_in_liquid( Atmosphere *A, const AtmosphereParameters Ap )
{
/* mass of volatiles in liquid */
PetscInt i;
PetscFunctionBeginUser;
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].mass_liquid = A->volatiles[i].x*A->Mliq;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_volatile_masses_in_solid( Atmosphere *A, const AtmosphereParameters Ap )
{
/* mass of volatiles in solid */
PetscInt i;
PetscFunctionBeginUser;
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].mass_solid = A->volatiles[i].x * Ap->volatile_parameters[i]->kdist * A->Msol;
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_volatile_masses_reactions( Atmosphere *A, const AtmosphereParameters Ap )
{
/* mass gain or loss due to reactions */
PetscInt i,j;
PetscScalar factor; /* normalisation, using the first volatile (reactant) in this chemical reaction */
PetscScalar massv;
PetscFunctionBeginUser;
/* initialise mass_reactions to zero for each volatile */
/* this is required, since these entries must be zero
without reactions to correctly compute volatile reservoirs */
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].mass_reaction = 0.0;
}
for (i=0; i<Ap->n_reactions; ++i) {
const PetscInt v0 = Ap->reaction_parameters[i]->volatiles[0];
/* by convention, first volatile is a reactant, so stoichiometry (hence factor) will be -ve */
factor = Ap->reaction_parameters[i]->stoichiometry[0] * Ap->volatile_parameters[v0]->molar_mass;
for (j=0; j<Ap->reaction_parameters[i]->n_volatiles; ++j) {
const PetscInt v = Ap->reaction_parameters[i]->volatiles[j];
massv = Ap->reaction_parameters[i]->stoichiometry[j] * Ap->volatile_parameters[v]->molar_mass;
massv /= factor; /* +ve for reactants, -ve for products */
A->volatiles[v].mass_reaction += massv * A->mass_reaction[i]; /* convention is mass loss of reactants, mass gain of products */
/* but regarding above, convention is arbitrary since mass_r[i] will switch sign to balance */
}
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_volume_mixing_ratios( Atmosphere *A, const AtmosphereParameters Ap )
{
PetscInt i;
PetscFunctionBeginUser;
A->molar_mass = 0.0;
for (i=0; i<Ap->n_volatiles; ++i) {
/* mixing ratio */
A->volatiles[i].mixing_ratio = A->volatiles[i].p / A->psurf;
/* atmosphere molar mass */
A->molar_mass += Ap->volatile_parameters[i]->molar_mass * A->volatiles[i].mixing_ratio;
}
PetscFunctionReturn(0);
}
PetscErrorCode set_reservoir_volatile_content( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, const ScalingConstants SC )
{
PetscErrorCode ierr;
PetscFunctionBeginUser;
/* if x is zero, the quantities below will also all
be set to zero, except the escape-related parameters that
are useful to compute even if the feedback of escape is not
actually included in the model */
/* order of these functions is very important! */
ierr = set_total_surface_pressure( A, Ap );CHKERRQ(ierr);
ierr = set_volatile_abundances_from_partial_pressure( A, Ap, SC );CHKERRQ(ierr);
ierr = set_volume_mixing_ratios( A, Ap );CHKERRQ(ierr);
ierr = set_volatile_masses_in_atmosphere( A, Ap );CHKERRQ(ierr);
ierr = set_volatile_masses_in_liquid( A, Ap );CHKERRQ(ierr);
ierr = set_volatile_masses_in_solid( A, Ap );CHKERRQ(ierr);
ierr = set_volatile_masses_reactions( A, Ap );CHKERRQ(ierr);
/* always set some escape terms, since then we can easily see the
magnitude of the effect, regardless of whether the feedback is
actually included for the volatile evolution */
ierr = set_jeans_escape( A, Ap, FC );CHKERRQ(ierr);
ierr = set_zahnle_escape( A, Ap );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_jeans_escape( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC )
{
PetscErrorCode ierr;
PetscInt i;
for (i=0; i<Ap->n_volatiles; ++i) {
ierr = set_column_density_volatile( A, Ap, FC, i );CHKERRQ(ierr);
ierr = set_jeans( A, Ap, FC, i );CHKERRQ(ierr);
ierr = set_Knudsen_number( A, Ap, i ); CHKERRQ(ierr);
ierr = set_R_thermal_escape( A, Ap, i ); CHKERRQ(ierr);
ierr = set_f_thermal_escape( A, i ); CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_f_zahnle_escape( Atmosphere *A, const AtmosphereParameters Ap, PetscInt i )
{
Volatile *V = &A->volatiles[i];
VolatileParameters const Vp = Ap->volatile_parameters[i];
PetscFunctionBeginUser;
/* Vp->R_zahnle_escape_value is set in parameters.c and contains all the
prefactors that are time-independent */
V->f_zahnle_escape = Vp->R_zahnle_escape_value * V->mixing_ratio;
PetscFunctionReturn(0);
}
static PetscErrorCode set_zahnle_escape( Atmosphere *A, const AtmosphereParameters Ap )
{
PetscErrorCode ierr;
PetscInt i;
for (i=0; i<Ap->n_volatiles; ++i) {
ierr = set_f_zahnle_escape( A, Ap, i );CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
PetscErrorCode set_atmosphere_emissivity_and_flux( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, const ScalingConstants SC )
{
PetscFunctionBeginUser;
switch( Ap->SURFACE_BC ){
case 1:
/* constant emissivity (greybody/blackbody) */
A->emissivity = Ap->emissivity0;
A->Fatm = get_grey_body_flux( A, Ap, FC );
break;
case 2:
/* Zahnle steam atmosphere parameterisation */
A->Fatm = get_steam_atmosphere_zahnle_1988_flux( A, SC );
A->emissivity = get_emissivity_from_flux( A, Ap, FC, A->Fatm );
break;
case 3:
/* Abe and Matsui (1985) */
A->emissivity = get_emissivity_abe_matsui( A, Ap );
A->Fatm = get_grey_body_flux( A, Ap, FC );
break;
case 4:
/* constant heat flux */
A->Fatm = Ap->surface_bc_value;
A->emissivity = get_emissivity_from_flux( A, Ap, FC, A->Fatm );
break;
case 5:
/* isothermal */
/* nothing to do, since bcs constrain the surface entropy
and entropy gradient directly */
break;
default:
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported SURFACE_BC value %d provided",Ap->SURFACE_BC);
break;
}
PetscFunctionReturn(0);
}
static PetscScalar get_grey_body_flux( const Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC )
{
PetscScalar Fsurf;
Fsurf = PetscPowScalar(A->tsurf,4.0)-PetscPowScalar(Ap->teqm,4.0);
Fsurf *= FC->STEFAN_BOLTZMANN * A->emissivity; /* Note emissivity may vary */
return Fsurf;
}
static PetscScalar get_steam_atmosphere_zahnle_1988_flux( const Atmosphere *A, const ScalingConstants SC )
{
PetscScalar Tsurf, Fsurf;
/* fit to Zahnle et al. (1988) from Solomatov and Stevenson (1993)
Eqn. 40. Expressed dimensionally so must convert here using
TEMP and FLUX scalings */
Tsurf = A->tsurf * SC->TEMP;
Fsurf = 1.5E2 + 1.02E-5 * PetscExpScalar(0.011*Tsurf);
Fsurf /= SC->FLUX; // non-dimensionalise
return Fsurf;
}
static PetscScalar get_emissivity_from_flux( const Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, PetscScalar flux )
{
PetscScalar emissivity;
emissivity = flux / FC->STEFAN_BOLTZMANN;
emissivity /= PetscPowScalar(A->tsurf,4.0)-PetscPowScalar(Ap->teqm,4.0);
return emissivity;
}
static PetscScalar get_emissivity_abe_matsui( Atmosphere *A, const AtmosphereParameters Ap )
{
PetscInt i;
PetscScalar emissivity;
A->tau = 0.0; // total surface optical depth
for (i=0; i<Ap->n_volatiles; ++i) {
Volatile *V = &A->volatiles[i];
/* optical depth at surface for this volatile */
V->tau = (3.0/2.0) * V->p / -(*Ap->gravity_ptr);
V->tau *= get_pressure_dependent_kabs( Ap, i );
/* total optical depth at surface */
A->tau += V->tau;
}
emissivity = 2.0 / (A->tau + 2.0);
return emissivity;
}
static PetscErrorCode JSON_add_atm_struct( Atmosphere *A, const AtmosphereParameters Ap, const FundamentalConstants FC, cJSON *data )
{
PetscErrorCode ierr;
PetscInt i;
PetscFunctionBeginUser;
/* only compute 1-D atmosphere structure for output */
ierr = set_atm_struct_tau( A );CHKERRQ(ierr);
ierr = set_atm_struct_temp( A, Ap, FC );CHKERRQ(ierr);
ierr = set_atm_struct_pressure( A, Ap );CHKERRQ(ierr);
ierr = set_atm_struct_depth( A, Ap, FC ); CHKERRQ(ierr);
/* write 1-D structure to JSON */
for (i=0;i<NUMATMSTRUCTVECS;++i){
cJSON *item;
DimensionalisableField curr = A->atm_struct[i];
ierr = DimensionalisableFieldToJSON(curr,&item);CHKERRQ(ierr);
cJSON_AddItemToObject(data,curr->name,item);
}
PetscFunctionReturn(0);
}
static PetscScalar get_pressure_dependent_kabs( const AtmosphereParameters Ap, PetscInt i )
{
/* absorption coefficient in the grey atmosphere is pressure-dependent
Abe and Matsui (1985), Eq. A21, A22, A23
*/
PetscScalar kabs;
kabs = PetscSqrtScalar( Ap->volatile_parameters[i]->kabs * -(*Ap->gravity_ptr) / (3.0*Ap->P0) );
return kabs;
}
PetscErrorCode JSON_add_atmosphere( DM dm, Parameters const P, Atmosphere *A, const char *name, cJSON *json )
{
PetscErrorCode ierr;
cJSON *data;
PetscScalar scaling, val;
FundamentalConstants const FC = P->fundamental_constants;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscInt v;
PetscFunctionBeginUser;
data = cJSON_CreateObject();
/* atmosphere properties only used by Abe and Matsui atmosphere model (1985) */
if(Ap->SURFACE_BC==3){
/* atmospheric structure */
ierr = JSON_add_atm_struct( A, Ap, FC, data );CHKERRQ(ierr);
/* optical depth, non-dimensional */
scaling = 1.0;
ierr = JSON_add_single_value_to_object(dm, scaling, "optical_depth", "None", A->tau, data);CHKERRQ(ierr);
}
if(Ap->SURFACE_BC!=5){
PetscScalar Teff, Tskin;
/* (effective) emissivity, non-dimensional */
scaling = 1.0;
ierr = JSON_add_single_value_to_object(dm, scaling, "emissivity", "None", A->emissivity, data);CHKERRQ(ierr);
/* net upward atmospheric flux */
scaling = SC->FLUX;
ierr = JSON_add_single_value_to_object(dm, scaling, "Fatm", "W m$^{-2}$", A->Fatm, data);CHKERRQ(ierr);
/* effective emitting temperature (T profile with scaled optical depth set to unity) */
Teff = A->Fatm / FC->STEFAN_BOLTZMANN + PetscPowScalar( Ap->teqm, 4.0 );
Teff = PetscPowScalar( Teff, 1.0/4.0 );
scaling = SC->TEMP;
ierr = JSON_add_single_value_to_object(dm, scaling, "Teff", "K", Teff, data);CHKERRQ(ierr);
Tskin = A->Fatm / (2*FC->STEFAN_BOLTZMANN) + PetscPowScalar( Ap->teqm, 4.0 );
Tskin = PetscPowScalar( Tskin, 1.0/4.0 );
ierr = JSON_add_single_value_to_object(dm, scaling, "Tskin", "K", Tskin, data);CHKERRQ(ierr);
}
/* total liquid mass of mantle, kg */
scaling = 4.0 * PETSC_PI * SC->MASS; // includes 4*PI for spherical geometry
ierr = JSON_add_single_value_to_object(dm, scaling, "mass_liquid", "kg", A->Mliq, data);CHKERRQ(ierr);
/* total solid mass of mantle, kg */
ierr = JSON_add_single_value_to_object(dm, scaling, "mass_solid", "kg", A->Msol, data);CHKERRQ(ierr);
/* total mass of mantle, kg (for sanity check) */
ierr = JSON_add_single_value_to_object(dm, scaling, "mass_mantle", "kg", *Ap->mantle_mass_ptr, data);CHKERRQ(ierr);
/* total mass of core, kg */
ierr = JSON_add_single_value_to_object(dm, scaling, "mass_core", "kg", P->coremass, data);CHKERRQ(ierr);
/* kg, without 4*pi */
val = SC->MASS * A->molar_mass * 1.0E3;
ierr = JSON_add_single_value_to_object(dm, 1.0, "molar_mass", "g/mol", val, data);CHKERRQ(ierr);
/* surface temperature, K */
scaling = SC->TEMP;
ierr = JSON_add_single_value_to_object(dm, scaling, "temperature_surface", "K", A->tsurf, data);CHKERRQ(ierr);
/* surface pressure, bar */
scaling = SC->PRESSURE / 1.0E5; /* bar */
ierr = JSON_add_single_value_to_object(dm, scaling, "pressure_surface", "bar", A->psurf, data);CHKERRQ(ierr);
/* oxygen fugacity */
if(Ap->OXYGEN_FUGACITY){
/* non-dimensional oxygen fugacity is pO2/psurf = volume mixing ratio of O2 */
scaling = 1.0;
ierr = JSON_add_single_value_to_object(dm, scaling, "fO2", "bar", A->fO2, data);CHKERRQ(ierr);
}
/* Volatiles */
for (v=0; v<Ap->n_volatiles; ++v) {
ierr = JSON_add_volatile(dm, P, Ap->volatile_parameters[v], &A->volatiles[v], A, Ap->volatile_parameters[v]->prefix, data ); CHKERRQ(ierr);
}
cJSON_AddItemToObject(json,name,data);
PetscFunctionReturn(0);
}
static PetscErrorCode JSON_add_volatile( DM dm, Parameters const P, VolatileParameters const VP, Volatile const *V, Atmosphere const *A, char const *name, cJSON *json )
{
PetscErrorCode ierr;
cJSON *data;
PetscScalar scaling, val;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscFunctionBeginUser;
data = cJSON_CreateObject();
/* parts-per-million (ppmw) */
scaling = SC->VOLATILE * 1.0E6; // VOLATILE to mass fraction, 1.0E6 to ppm
/* initial volatile (ppmw) */
ierr = JSON_add_single_value_to_object(dm, scaling, "initial_ppmw", "ppmw", VP->initial_total_abundance, data);CHKERRQ(ierr);
/* volatile in liquid mantle (ppmw) */
ierr = JSON_add_single_value_to_object(dm, scaling, "liquid_ppmw", "ppmw", V->x, data);CHKERRQ(ierr);
/* volatile in solid mantle (ppmw) */
ierr = JSON_add_single_value_to_object(dm, scaling, "solid_ppmw", "ppmw", V->x*VP->kdist, data);CHKERRQ(ierr);
/* ocean moles */
scaling = 1.0;
ierr = JSON_add_single_value_to_object(dm, scaling, "initial_ocean_moles", "ocean moles", VP->initial_ocean_moles, data);CHKERRQ(ierr);
/* kilograms (kg) */
scaling = SC->VOLATILE * 4.0 * PETSC_PI * SC->MASS;
/* initial volatile (kg) */
ierr = JSON_add_single_value_to_object(dm, scaling, "initial_kg", "kg", VP->initial_total_abundance*(*Ap->mantle_mass_ptr), data);CHKERRQ(ierr);
/* volatile in liquid mantle (kg) */
ierr = JSON_add_single_value_to_object(dm, scaling, "liquid_kg", "kg", V->mass_liquid, data);CHKERRQ(ierr);
/* volatile in solid mantle (kg) */
ierr = JSON_add_single_value_to_object(dm, scaling, "solid_kg", "kg", V->mass_solid, data);CHKERRQ(ierr);
/* volatile in atmosphere (kg) */
ierr = JSON_add_single_value_to_object(dm, scaling, "atmosphere_kg", "kg", V->mass_atmos, data);CHKERRQ(ierr);
/* volatile exchanges due to reactions (kg) */
ierr = JSON_add_single_value_to_object(dm, scaling, "reaction_kg", "kg", V->mass_reaction, data);CHKERRQ(ierr);
/* physical reservoir size */
val = V->mass_liquid + V->mass_solid + V->mass_atmos;
ierr = JSON_add_single_value_to_object(dm, scaling, "physical_kg", "kg", val, data);CHKERRQ(ierr);
/* grams (g), without 4*pi */
val = SC->MASS * VP->molar_mass * 1.0E3;
ierr = JSON_add_single_value_to_object(dm, 1.0, "molar_mass", "g/mol", val, data);CHKERRQ(ierr);
/* bar (bar) */
/* volatile in atmosphere (bar) */
scaling = SC->PRESSURE / 1.0E5; /* bar */
ierr = JSON_add_single_value_to_object(dm, scaling, "atmosphere_bar", "bar", V->p, data);CHKERRQ(ierr);
scaling = SC->PRESSURE / 1.0E5 / SC->TIMEYRS; /* bar / year */
ierr = JSON_add_single_value_to_object(dm, scaling, "dpdt", "bar/yr", V->dpdt, data);CHKERRQ(ierr);
/* area */
scaling = 1.0 / (SC->AREA * 1.0E4); // 1/cm^2
ierr = JSON_add_single_value_to_object(dm, scaling, "column_density", "1/cm^2", V->column_density, data);CHKERRQ(ierr);
/* non-dimensional */
scaling = 1.0;
ierr = JSON_add_single_value_to_object(dm, scaling, "optical_depth", "None", V->tau, data);CHKERRQ(ierr);
ierr = JSON_add_single_value_to_object(dm, scaling, "mixing_ratio", "None", V->mixing_ratio, data);CHKERRQ(ierr);
ierr = JSON_add_single_value_to_object(dm, scaling, "jeans", "None", V->jeans, data);CHKERRQ(ierr);
ierr = JSON_add_single_value_to_object(dm, scaling, "Knudsen", "None", V->Knudsen, data);CHKERRQ(ierr);
ierr = JSON_add_single_value_to_object(dm, scaling, "R_thermal_escape", "None", V->R_thermal_escape, data);CHKERRQ(ierr);
ierr = JSON_add_single_value_to_object(dm, scaling, "f_thermal_escape", "None", V->f_thermal_escape, data);CHKERRQ(ierr);
cJSON_AddItemToObject(json,name,data);
PetscFunctionReturn(0);
}
PetscErrorCode objective_function_volatile_evolution( SNES snes, Vec x, Vec f, void *ptr)
{
PetscErrorCode ierr;
const PetscScalar *xx;
PetscScalar *ff;
Ctx *E = (Ctx*) ptr;
PetscInt i;
const PetscScalar *dmrdt;
Atmosphere *A = &E->atmosphere;
Parameters const P = E->parameters;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscFunctionBeginUser;
ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].dpdt = xx[i];
}
dmrdt = &xx[Ap->n_volatiles];
ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
for (i=0; i<Ap->n_volatiles; ++i) {
ff[i] = get_dpdt( A, Ap, i, dmrdt, SC );
}
/* chemical equilibrium constraints */
for (i=0; i<Ap->n_reactions; ++i) {
PetscScalar dQpdt, dQrdt, Qr, dGdt, G, log10G, dlog10GdT;
Qr = get_reaction_quotient_reactants( Ap->reaction_parameters[i], A );
//Qp = get_reaction_quotient_products( Ap->reaction_parameters[i], A );
dQpdt = get_reaction_quotient_products_time_derivative( Ap->reaction_parameters[i], A, Ap );
dQrdt = get_reaction_quotient_reactants_time_derivative( Ap->reaction_parameters[i], A, Ap );
log10G = get_log10_modified_equilibrium_constant( Ap->reaction_parameters[i], A->tsurf, A, SC );
dlog10GdT = get_dlog10GdT( Ap->reaction_parameters[i], A->tsurf, A );
/* deriv of Qp/Qr - G */
//ff[Ap->n_volatiles + i] = -Qp/PetscPowScalar(Qr,2.0) * dQrdt;
//ff[Ap->n_volatiles + i] += 1.0/Qr * dQpdt;
/* instead, multiply by Qr to avoid division. This simplifies
the evaluation of the derivative, and is fine if p's are
scaled to around unity */
ff[Ap->n_volatiles + i] = dQpdt;
/* Modified equilibrium constant */
G = PetscPowScalar( 10.0, log10G );
/* dG/dlog10G * dlog10G/dT * dT/dt */
dGdt = G * PetscLogReal( 10.0 ) * dlog10GdT * A->dtsurfdt;
ff[Ap->n_volatiles + i] -= dGdt*Qr + G*dQrdt;
}
ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
PetscScalar get_dpdt( Atmosphere *A, const AtmosphereParameters Ap, PetscInt i, const PetscScalar *dmrdt, const ScalingConstants SC )
{
PetscErrorCode ierr;
PetscScalar out, out2, massv, f_thermal_escape, f_constant_escape, f_zahnle_escape;
PetscInt j,k;
VolatileParameters const Vp = Ap->volatile_parameters[i];
Volatile *V = &A->volatiles[i];
PetscScalar pbar=0, Tkel=0;
/* remember that to this point, V->f_thermal_escape is always
computed so we can see the (potential) role of thermal escape,
but it is not necessarily used in the calculation */
if(Ap->THERMAL_ESCAPE){
f_thermal_escape = V->f_thermal_escape;
}
else{
/* thermal escape of unity is necessary to scale the source (outgassing) rate,
but the other modifiers due to Jean's escape are not included */
f_thermal_escape = 1.0;
}
/* constant escape, non-thermal (Jean's) contribution */
if(Ap->CONSTANT_ESCAPE){
f_constant_escape = Vp->constant_escape_value;
}
else{
f_constant_escape = 0.0;
}
/* Zahnle et al. (2019, Eqn 3, Fig 5), escape */
if(Ap->ZAHNLE_ESCAPE){
f_zahnle_escape = V->f_zahnle_escape;
}
else{
f_zahnle_escape = 0.0;
}
out2 = 0.0;
ierr = set_total_surface_pressure_time_derivative( A, Ap );CHKERRQ(ierr);
/* next part concerns d/dt relating to the mean molar mass of the atmosphere */
for (j=0; j<Ap->n_volatiles; ++j) {
out = -A->dpsurfdt * A->volatiles[j].p / A->psurf;
out += A->volatiles[j].dpdt;
out *= Ap->volatile_parameters[j]->molar_mass;
out2 += out;
}
out2 *= -V->p / (A->psurf * PetscSqr(A->molar_mass));
/* second part of atmosphere derivative */
out2 += ( 1.0 / A->molar_mass ) * V->dpdt;
/* multiply by prefactors */
out2 *= (1.0 / (*Ap->VOLATILE_ptr)) * PetscSqr(*Ap->radius_ptr) * Vp->molar_mass / -(*Ap->gravity_ptr); // note negative gravity
/* thermal (Jean's) escape correction */
out2 *= f_thermal_escape;
/* non-thermal escape contribution */
out2 += f_constant_escape;
/* Zahnle et al. (2019) escape contribution */
out2 += f_zahnle_escape;
/* chemical reactions. Loop through all reactions, check if they involve
this volatile, and if so, add a term */
for (j=0; j<Ap->n_reactions; ++j) {
/* normalisation, using the first volatile (reactant) in this chemical reaction */
const PetscInt v0 = Ap->reaction_parameters[j]->volatiles[0];
for (k=0; k<Ap->reaction_parameters[j]->n_volatiles; ++k) {
const PetscInt v = Ap->reaction_parameters[j]->volatiles[k];
/* i is the current volatile of interest */
if (v==i) {
massv = dmrdt[j];
massv *= Ap->reaction_parameters[j]->stoichiometry[k] / Ap->reaction_parameters[j]->stoichiometry[0];
massv *= Ap->volatile_parameters[v]->molar_mass / Ap->volatile_parameters[v0]->molar_mass;
out2 += massv;
}
}
}
switch( Vp->SOLUBILITY ){