-
Notifications
You must be signed in to change notification settings - Fork 0
/
ic.c
1024 lines (822 loc) · 38.9 KB
/
ic.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 "atmosphere.h"
#include "bc.h"
#include "energy.h"
#include "eos.h"
#include "ic.h"
#include "matprop.h"
#include "monitor.h"
#include "parameters.h"
#include "reaction.h"
#include "util.h"
/* interior ic */
static PetscErrorCode set_ic_interior( Ctx *, Vec );
static PetscErrorCode set_ic_interior_default( Ctx *, Vec );
static PetscErrorCode set_ic_interior_entropy( Ctx *, Vec );
static PetscErrorCode set_ic_interior_from_file( Ctx *, Vec );
static PetscErrorCode set_ic_interior_from_phase_boundary( Ctx *, Vec );
static PetscErrorCode set_start_time_from_file( Parameters , const char * );
static PetscErrorCode set_start_stepmacro_from_file( Parameters , const char * );
/* atmosphere ic */
static PetscErrorCode set_ic_atmosphere( Ctx *, Vec );
static PetscErrorCode set_ic_atmosphere_pseudo_volatiles( Ctx *, Vec );
static PetscErrorCode set_ic_atmosphere_from_ocean_moles( Ctx *E, Vec sol );
static PetscErrorCode set_ic_atmosphere_from_initial_total_abundance( Ctx *E, Vec sol );
static PetscErrorCode set_ic_atmosphere_from_file( Ctx *, Vec );
static PetscErrorCode set_ic_atmosphere_from_partial_pressure( Ctx *, Vec );
static PetscErrorCode solve_for_initial_partial_pressure( Ctx * );
/* general */
static PetscErrorCode set_ic_from_file( Ctx *, Vec, const char *, const PetscInt *, PetscInt );
static PetscErrorCode objective_function_initial_partial_pressure( SNES, Vec, Vec, void *);
static PetscErrorCode conform_atmosphere_parameters_to_ic( Ctx * );
static PetscErrorCode print_ocean_masses( Ctx * );
/* main function to set initial condition */
PetscErrorCode set_initial_condition( Ctx *E, Vec sol)
{
PetscErrorCode ierr;
Parameters const P = E->parameters;
Atmosphere *A = &E->atmosphere;
AtmosphereParameters Ap = P->atmosphere_parameters;
FundamentalConstants FC = P->fundamental_constants;
ScalingConstants SC = P->scaling_constants;
PetscFunctionBeginUser;
/* interior initial condition */
ierr = set_ic_interior( E, sol ); CHKERRQ(ierr);
/* need some interface quantities for setting atmosphere ic */
ierr = set_interior_atmosphere_interface_from_surface_entropy( E ); CHKERRQ(ierr);
/* atmosphere initial condition */
/* must be after interior initial condition */
ierr = set_ic_atmosphere( E, sol ); CHKERRQ(ierr);
/* for surface bc, need an estimate of A->Fatm */
ierr = set_atmosphere_emissivity_and_flux( A, Ap, FC, SC );
/* P->t0 is set in parameters.c */
/* time is needed for the radiogenic energy input */
/* if Ap->SURFACE_BC_ACC is set, then the surface radiation is solved
for every time step. But sometimes we just want to solve the balance
once for the IC, so define Ap->SURFACE_BC_ACC_IC for this purpose */
if( (P->ic_steady_state_energy) || (Ap->SURFACE_BC_ACC)){
ierr = solve_for_surface_radiation_balance( E, P->t0 );CHKERRQ(ierr);
/* we also solve for the steady-state initially, by ensuring the energy
flows in the interior are pretty much equal to the surface energy
flow. This should help the time-stepper to get going a bit faster
and avoid some clunky spurious values from applying an unphysical
initial condition */
ierr = solve_for_steady_state_energy_interior( E, P->t0 );CHKERRQ(ierr);
ierr = set_solution_from_entropy( E, sol );CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
/* initial condition of interior */
static PetscErrorCode set_ic_interior( Ctx *E, Vec sol)
{
PetscErrorCode ierr;
Parameters const P = E->parameters;
PetscFunctionBeginUser;
if (P->IC_INTERIOR==1){
ierr = set_ic_interior_default( E, sol ); CHKERRQ(ierr);
}
else if (P->IC_INTERIOR==2){
ierr = set_ic_interior_from_file( E, sol ); CHKERRQ(ierr);
ierr = set_start_time_from_file( P, P->ic_interior_filename );CHKERRQ(ierr);
ierr = set_start_stepmacro_from_file( P, P->ic_interior_filename );CHKERRQ(ierr);
}
else if (P->IC_INTERIOR==3){
ierr = set_ic_interior_from_phase_boundary( E, sol ); CHKERRQ(ierr);
}
else{
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported IC_INTERIOR value %d provided",P->IC_INTERIOR);
}
/* this copies the sol Vecs to E */
/* the reason we set sol above, rather than Vecs in struct, is
because the function below is called by the RHS and must take
sol as the input. This now sets the Vecs in E */
ierr = set_entropy_from_solution(E, sol); CHKERRQ(ierr);
/* option to set initial core-mantle boundary entropy */
if( P->ic_core_entropy > 0.0 ){
ierr = set_cmb_entropy_constant( E );CHKERRQ(ierr);
}
/* option to set initial surface entropy */
if( P->ic_surface_entropy > 0.0 ){
ierr = set_surface_entropy_constant( E );CHKERRQ(ierr);
}
/* Now the Vecs in E are set and consistent, clone them
to the sol Vec which is what is actually used by the
time-stepper */
ierr = set_solution_from_entropy(E, sol);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_interior_default( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
PetscFunctionBeginUser;
ierr = set_ic_interior_entropy( E, sol ); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_interior_from_file( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
Parameters const P = E->parameters;
/* subdomains to use, i.e. dS/dxi and S0 */
PetscInt const arr[2] = {0, 1};
PetscFunctionBeginUser;
ierr = set_ic_from_file( E, sol, P->ic_interior_filename, arr, 2 ); CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_atmosphere_from_file( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
Parameters const P = E->parameters;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscInt const arr1[1] = {2};
PetscInt const arr2[1] = {3};
PetscFunctionBeginUser;
if(Ap->n_volatiles){
ierr = set_ic_from_file( E, sol, Ap->ic_atmosphere_filename, arr1, 1 ); CHKERRQ(ierr);
}
if(Ap->n_reactions){
ierr = set_ic_from_file( E, sol, Ap->ic_atmosphere_filename, arr2, 1 ); CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_interior_entropy( Ctx *E, Vec sol )
{
/* set initial entropy gradient to constant for all basic nodes */
PetscErrorCode ierr;
PetscInt i;
PetscScalar S0;
Parameters const P = E->parameters;
Mesh const *M = &E->mesh;
Vec dSdxi_b;
Vec *subVecs;
PetscFunctionBeginUser;
ierr = PetscMalloc1(E->numFields,&subVecs);CHKERRQ(ierr);
ierr = DMCompositeGetAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
dSdxi_b = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_DSDXI_B]];
/* set entropy gradient and map to mass coordinates */
ierr = VecSet(dSdxi_b, P->ic_dsdr );CHKERRQ(ierr);
ierr = VecPointwiseDivide(dSdxi_b,dSdxi_b,M->dxidr_b);
/* code to test surface bc (not required other than for debugging) */
//PetscInt const ind0=0,ind1=1;
//PetscScalar val1;
//ierr = VecGetValues(dSdxi_b,1,&ind1,&val1);CHKERRQ(ierr);
//ierr = VecSetValue(dSdxi_b,ind0,val1,INSERT_VALUES);CHKERRQ(ierr);
/* set entropy at top of adiabat */
S0 = P->ic_adiabat_entropy;
ierr = VecSetValue(subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_S0]],0,S0,INSERT_VALUES);CHKERRQ(ierr);
for (i=0; i<E->numFields; ++i) {
ierr = VecAssemblyBegin(subVecs[i]);CHKERRQ(ierr);
ierr = VecAssemblyEnd(subVecs[i]);CHKERRQ(ierr);
}
ierr = DMCompositeRestoreAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
ierr = PetscFree(subVecs);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_start_time_from_file( Parameters P , const char * filename )
{
PetscErrorCode ierr;
cJSON *json=NULL, *time;
PetscFunctionBeginUser;
ierr = read_JSON_file_to_JSON_object( filename, &json );
/* time from this restart JSON must be passed to the time stepper
to continue the integration and keep track of absolute time */
time = cJSON_GetObjectItem(json,"time");
P->t0 = time->valuedouble;
cJSON_Delete( json );
PetscFunctionReturn(0);
}
static PetscErrorCode set_start_stepmacro_from_file( Parameters P, const char * filename )
{
PetscErrorCode ierr;
cJSON *json=NULL, *stepmacro;
PetscFunctionBeginUser;
ierr = read_JSON_file_to_JSON_object( filename, &json );
/* stepmacro from this restart JSON must be passed to the time stepper
to continue the integration */
stepmacro = cJSON_GetObjectItem(json,"step");
P->stepmacro = stepmacro->valueint;
cJSON_Delete( json );
PetscFunctionReturn(0);
}
PetscErrorCode read_JSON_file_to_JSON_object( const char * filename, cJSON ** json )
{
FILE *fp;
long length;
char *buffer = 0;
PetscFunctionBeginUser;
fp = fopen( filename, "r" );
if(fp==NULL) {
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Could not open file %s",filename);
}
/* read file to zero terminated string */
/* https://stackoverflow.com/questions/174531/how-to-read-the-content-of-a-file-to-a-string-in-c */
if (fp){
fseek(fp,0,SEEK_END);
length = ftell(fp);
fseek(fp,0,SEEK_SET);
buffer = malloc(length+1);
if (buffer){
if(fread(buffer,1,length,fp) != length) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_FILE_READ,"fread() error");
} else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_MEM,"malloc error");
fclose(fp);
buffer[length] = '\0';
}
/* for future reference for returning pointers */
/* https://stackoverflow.com/questions/2520777/pointer-as-second-argument-instead-of-returning-pointer */
*json = cJSON_Parse( buffer );
free( buffer );
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_from_file( Ctx *E, Vec sol, const char * filename, const PetscInt *arr, PetscInt size_arr )
{
/* reads an initial condition from a previously output JSON file
to enable restarting */
PetscErrorCode ierr;
cJSON *json=NULL, *solution, *subdomain, *values, *data, *item;
char *item_str;
PetscInt i, j, subdomain_num;
PetscScalar val = 0;
Vec invec, *subVecs;
#if (defined PETSC_USE_REAL___FLOAT128)
char val_str[30];
#endif
PetscFunctionBeginUser;
ierr = PetscMalloc1(E->numFields,&subVecs);CHKERRQ(ierr);
ierr = DMCompositeGetAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
ierr = read_JSON_file_to_JSON_object( filename, &json );
solution = cJSON_GetObjectItem(json,"solution");
subdomain = cJSON_GetObjectItem(solution,"subdomain data");
/* loop over all subdomains */
cJSON_ArrayForEach( data, subdomain )
{
subdomain = cJSON_GetObjectItem( data, "subdomain" );
subdomain_num = subdomain->valueint;
values = cJSON_GetObjectItem( data, "values" );
/* could break if ordering of subdomains changes */
if (subdomain_num == 0){
invec = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_DSDXI_B]];
}
else if (subdomain_num == 1){
invec = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_S0]];
}
else if (subdomain_num == 2){
invec = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_MO_VOLATILES]];
}
else if (subdomain_num == 3){
invec = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_MO_REACTIONS]];
}
else {
SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"unexpected number of subdomains");
}
/* determine if this subdomain is required for the IC by checking
against the listed subdomains in arr */
for( j=0; j<size_arr; j++ ){
if(subdomain_num == arr[j]){
for (i=0; i<cJSON_GetArraySize(values); i++ ){
item = cJSON_GetArrayItem( values, i );
item_str = item->valuestring;
#if (defined PETSC_USE_REAL___FLOAT128)
sscanf( item_str, "%s", val_str );
val = strtoflt128(val_str, NULL);
#else
sscanf( item_str, "%lf", &val );
#endif
/* add value to vec */
VecSetValue( invec, i, val, INSERT_VALUES );CHKERRQ(ierr);
}
}
}
}
/* all vec assembly is done here */
for (i=0; i<E->numFields; ++i) {
ierr = VecAssemblyBegin(subVecs[i]);CHKERRQ(ierr);
ierr = VecAssemblyEnd(subVecs[i]);CHKERRQ(ierr);
}
ierr = DMCompositeRestoreAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
ierr = PetscFree(subVecs);CHKERRQ(ierr);
cJSON_Delete( json );
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_interior_from_phase_boundary( Ctx *E, Vec sol )
{
/* entropy tracks a phase boundary (usually the solidus) below a cutoff value, and is
equal to the cutoff value above */
PetscErrorCode ierr;
Parameters const P = E->parameters;
Mesh const *M = &E->mesh;
PetscScalar S0,S1,S2;
const PetscScalar *arr_pres_s;
PetscScalar *arr_dSdxi_b, *arr_xi_s;
PetscInt i, ihi_b, ilo_b, w_b;
DM da_s = E->da_s, da_b=E->da_b;
Vec dSdxi_b, *subVecs;
PetscFunctionBeginUser;
ierr = PetscMalloc1(E->numFields,&subVecs);CHKERRQ(ierr);
ierr = DMCompositeGetAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
dSdxi_b = subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_DSDXI_B]];
/* for looping over basic nodes */
ierr = DMDAGetCorners(da_b,&ilo_b,0,0,&w_b,0,0);CHKERRQ(ierr);
ihi_b = ilo_b + w_b;
ierr = DMDAVecGetArray(da_b,dSdxi_b,&arr_dSdxi_b);CHKERRQ(ierr);
ierr = DMDAVecGetArrayRead(da_s,M->xi_s,&arr_xi_s);CHKERRQ(ierr);
ierr = DMDAVecGetArrayRead(E->da_s,M->pressure_s,&arr_pres_s);CHKERRQ(ierr);
/* first and last points not updated by loop below */
arr_dSdxi_b[0] = 0.0;
arr_dSdxi_b[ihi_b-1] = 0.0;
for(i=ilo_b+1; i<ihi_b-1; ++i){
/* assumes the relevant phase boundary is in slot 0, which it will be for a single
phase system (although still potential for a bug here if you want to initialise a
two phase system from a phase boundary ic) */
/* get staggered nodes values either side of the basic node */
ierr = EOSGetPhaseBoundary( E->parameters->eos_phases[0], arr_pres_s[i], &S1, NULL);CHKERRQ(ierr);
S1 = PetscMin( S1, P->ic_adiabat_entropy );
ierr = EOSGetPhaseBoundary( E->parameters->eos_phases[0], arr_pres_s[i-1], &S2, NULL);CHKERRQ(ierr);
S2 = PetscMin( S2, P->ic_adiabat_entropy );
/* now compute gradient at basic node */
arr_dSdxi_b[i] = S1 - S2;
arr_dSdxi_b[i] /= arr_xi_s[i] - arr_xi_s[i-1];
}
/* for last point, assume same gradient as basic node above */
arr_dSdxi_b[ihi_b-1] = arr_dSdxi_b[ihi_b-2];
/* set entropy at top staggered node */
ierr = EOSGetPhaseBoundary( E->parameters->eos_phases[0], arr_pres_s[0], &S0, NULL);CHKERRQ(ierr);
S0 = PetscMin( S0, P->ic_adiabat_entropy );
ierr = VecSetValue(subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_S0]],0,S0,INSERT_VALUES);CHKERRQ(ierr);
ierr = VecAssemblyBegin(subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_S0]]);CHKERRQ(ierr);
ierr = VecAssemblyEnd(subVecs[E->solutionSlots[SPIDER_SOLUTION_FIELD_S0]]);CHKERRQ(ierr);
/* restore basic node gradient array and sol Vec */
ierr = DMDAVecRestoreArray(da_b,dSdxi_b,&arr_dSdxi_b);CHKERRQ(ierr);
ierr = DMCompositeRestoreAccessArray(E->dm_sol,sol,E->numFields,NULL,subVecs);CHKERRQ(ierr);
ierr = PetscFree(subVecs);CHKERRQ(ierr);
ierr = DMDAVecRestoreArrayRead(da_s,M->xi_s,&arr_xi_s);CHKERRQ(ierr);
ierr = DMDAVecRestoreArrayRead(da_s,M->pressure_s,&arr_pres_s);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
/* initial condition of atmosphere */
static PetscErrorCode set_ic_atmosphere( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
Parameters const P = E->parameters;
FundamentalConstants const FC = P->fundamental_constants;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
Atmosphere *A = &E->atmosphere;
PetscFunctionBeginUser;
if(Ap->n_volatiles){
/* for pseudo-volatiles we get the atmospheric partial pressures from
lookup, based on the surface temperature */
if(Ap->PSEUDO_VOLATILES){
ierr = set_ic_atmosphere_pseudo_volatiles( E, sol );CHKERRQ(ierr);
}
/* set atmosphere IC from initial total abundance */
else if(Ap->IC_ATMOSPHERE==1){
ierr = set_ic_atmosphere_from_initial_total_abundance( E, sol );CHKERRQ(ierr);
}
/* read atmosphere IC from file (could be a different file to
the interior IC) */
else if (Ap->IC_ATMOSPHERE==2){
/* this only sets the sol Vec, and not the A->volatiles[i].x */
ierr = set_ic_atmosphere_from_file( E, sol );CHKERRQ(ierr);
}
/* set IC from partial pressure in the atmosphere */
else if (Ap->IC_ATMOSPHERE==3){
ierr = set_ic_atmosphere_from_partial_pressure( E, sol );CHKERRQ(ierr);
}
/* set IC from number of ocean moles */
else if (Ap->IC_ATMOSPHERE==4){
ierr = set_ic_atmosphere_from_ocean_moles( E, sol );CHKERRQ(ierr);
}
else{
SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported IC_ATMOSPHERE value %d provided",Ap->IC_ATMOSPHERE);
}
/* all of the above are guaranteed to set Vec sol, for the
partial pressures and the mass reaction terms. Now also
ensure that A->volatiles[i].p and A->mass_reaction[i]
conform to the Vec sol, since these are subsequently
used to "correct" the mass reaction to zero */
ierr = set_partial_pressures_from_solution( E, sol );CHKERRQ(ierr);
/* ensure all atmosphere quantities are consistent with current
solution */
ierr = set_reservoir_volatile_content( A, Ap, FC, SC );CHKERRQ(ierr);
/* might mess up pseudo-volatiles? */
ierr = conform_atmosphere_parameters_to_ic( E );CHKERRQ(ierr);
ierr = set_solution_from_partial_pressures( E, sol );CHKERRQ(ierr);
}
/* with Ap->n_volatiles=0 there are no entries in
the solution vector to initialise to zero */
PetscFunctionReturn(0);
}
static PetscErrorCode conform_atmosphere_parameters_to_ic( Ctx *E )
{
/* this function has authorisation to update the parameters struct
to ensure that parameters obey constraints imposed by the
initial condition */
PetscErrorCode ierr;
PetscInt i;
PetscScalar mass;
Parameters P = E->parameters;
Atmosphere *A = &E->atmosphere;
AtmosphereParameters Ap = P->atmosphere_parameters;
FundamentalConstants const FC = P->fundamental_constants;
ScalingConstants const SC = P->scaling_constants;
PetscFunctionBeginUser;
/* prior to this function, A->volatiles[i].p and A->mass_reaction[i]
are updated */
for (i=0; i<Ap->n_volatiles; ++i){
/* to this point, we have solved for A->volatiles[i].p.
we must now ensure that other parameters are self-consistent
with this pressure. This might over-write some parameters
that were previously used to determine p. */
/* even if we restart, we can conform these parameters */
mass = A->volatiles[i].mass_liquid + A->volatiles[i].mass_solid + A->volatiles[i].mass_atmos + A->volatiles[i].mass_reaction;
Ap->volatile_parameters[i]->initial_total_abundance = mass / (*Ap->mantle_mass_ptr);
/* initial partial pressure */
Ap->volatile_parameters[i]->initial_atmos_pressure = A->volatiles[i].p;
/* initial_ocean_moles */
Ap->volatile_parameters[i]->initial_ocean_moles = Ap->volatile_parameters[i]->initial_total_abundance;
Ap->volatile_parameters[i]->initial_ocean_moles *= (*Ap->mantle_mass_ptr) * SC->VOLATILE * 4.0 * PETSC_PI;
Ap->volatile_parameters[i]->initial_ocean_moles /= Ap->volatile_parameters[i]->molar_mass;
Ap->volatile_parameters[i]->initial_ocean_moles /= FC->OCEAN_MOLES;
}
/* if not restarting */
if( Ap->IC_ATMOSPHERE != 2 ){
/* since we updated the parameters above to adhere to the equilibrium
state, we reset the mass_reaction to zero since our initial
abundances are referenced to equilibrium at the initial
interior and surface conditions (e.g., A->tsurf) */
for(i=0; i<Ap->n_reactions; ++i){
A->mass_reaction[i] = 0.0;
}
}
ierr = print_ocean_masses( E );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode print_ocean_masses( Ctx *E )
{
/* useful to report the equivalent present-day ocean mass
1.4E21 kg of H2O (Olson and Sharp, 2018)
1.55E20 kg of H2 (Olson and Sharp, 2018) */
/* this is not general, and will break if other species contain H2,
e.g., CH4. Ideally, we should loop over all volatiles and identify
those that include H2 */
PetscErrorCode ierr;
PetscInt i;
/* for an arbitrary number of volatile species, would be better to somehow
automate this output, rather than making it specific to reduced and oxidised
phases of hydrogen and oxygen */
PetscBool FLAG_H2 = PETSC_FALSE;
PetscBool FLAG_H2O = PETSC_FALSE;
PetscBool FLAG_CO = PETSC_FALSE;
PetscBool FLAG_CO2 = PETSC_FALSE;
PetscBool FLAG_CH4 = PETSC_FALSE;
PetscScalar mass_H2 = 0.0, molar_mass_H2 = 0.0, tmass_H2 = 0.0, p_H2 = 0.0;
PetscScalar mass_H2O = 0.0, molar_mass_H2O = 0.0, tmass_H2O = 0.0, p_H2O = 0.0;
PetscScalar mass_CO = 0.0, molar_mass_CO = 0.0, tmass_CO = 0.0, p_CO = 0.0;
PetscScalar mass_CO2 = 0.0, molar_mass_CO2 = 0.0, tmass_CO2 = 0.0, p_CO2 = 0.0;
PetscScalar mass_CH4 = 0.0, molar_mass_CH4 = 0.0; //, tmass_CH4 = 0.0, p_CH4 = 0.0;
Parameters P = E->parameters;
AtmosphereParameters Ap = P->atmosphere_parameters;
ScalingConstants const SC = P->scaling_constants;
PetscScalar scaling = SC->VOLATILE * 4.0 * PETSC_PI * SC->MASS;
PetscScalar scaling2 = (1.0/SC->VOLATILE) * PetscSqr(*Ap->radius_ptr) / -(*Ap->gravity_ptr);
PetscFunctionBeginUser;
/* this is ugly, but otherwise the flags get reset if they appear within
the same loop over volatiles */
for(i=0; i<Ap->n_volatiles; ++i){
PetscStrcmp( Ap->volatile_parameters[i]->prefix, "H2", &FLAG_H2 );
if ( FLAG_H2 ){
mass_H2 = Ap->volatile_parameters[i]->initial_total_abundance;
molar_mass_H2 = Ap->volatile_parameters[i]->molar_mass;
break;
}
}
for(i=0; i<Ap->n_volatiles; ++i){
PetscStrcmp( Ap->volatile_parameters[i]->prefix, "H2O", &FLAG_H2O );
if ( FLAG_H2O ){
mass_H2O = Ap->volatile_parameters[i]->initial_total_abundance;
molar_mass_H2O = Ap->volatile_parameters[i]->molar_mass;
break;
}
}
for(i=0; i<Ap->n_volatiles; ++i){
PetscStrcmp( Ap->volatile_parameters[i]->prefix, "CO", &FLAG_CO );
if ( FLAG_CO ){
mass_CO = Ap->volatile_parameters[i]->initial_total_abundance;
molar_mass_CO = Ap->volatile_parameters[i]->molar_mass;
break;
}
}
for(i=0; i<Ap->n_volatiles; ++i){
PetscStrcmp( Ap->volatile_parameters[i]->prefix, "CO2", &FLAG_CO2 );
if ( FLAG_CO2 ){
mass_CO2 = Ap->volatile_parameters[i]->initial_total_abundance;
molar_mass_CO2 = Ap->volatile_parameters[i]->molar_mass;
break;
}
}
for(i=0; i<Ap->n_volatiles; ++i){
PetscStrcmp( Ap->volatile_parameters[i]->prefix, "CH4", &FLAG_CH4 );
if ( FLAG_CH4 ){
mass_CH4 = Ap->volatile_parameters[i]->initial_total_abundance;
molar_mass_CH4 = Ap->volatile_parameters[i]->molar_mass;
break;
}
}
if( FLAG_H2 || FLAG_H2O || FLAG_CO || FLAG_CO2 || FLAG_CH4 ){
ierr = PetscPrintf(PETSC_COMM_WORLD,"\n**************** Volatile content **************\n");CHKERRQ(ierr);
}
/* the output below should report the same value, but it is a useful
check that everything is working as expected */
/* H2 */
if( FLAG_H2 ){
tmass_H2 = mass_H2;
if ( FLAG_H2O ){
/* equivalent mass of H2 in H2O */
tmass_H2 += mass_H2O * (molar_mass_H2 / molar_mass_H2O );
}
if ( FLAG_CH4 ){
/* equivalent mass of H2 in CH4 */
tmass_H2 += mass_CH4 * (2 * molar_mass_H2 / molar_mass_CH4 );
}
tmass_H2 *= (*Ap->mantle_mass_ptr); /* total non-dimensional mass */
p_H2 = tmass_H2 / scaling2; /* equivalent surface pressure */
tmass_H2 *= scaling; /* total physical mass */
tmass_H2 /= 1.55E20; /* non-dimensionalise according to ocean mass of H2 */
p_H2 *= SC->PRESSURE / 1.0E5; /* to bar */
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent present-day mass of ocean water from H2 (non-dimensional)",(double)tmass_H2);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent atmospheric pressure of H2 (bar)",(double)p_H2);CHKERRQ(ierr);
}
/* H2O */
if( FLAG_H2O ){
tmass_H2O = mass_H2O;
if ( FLAG_H2 ){
/* equivalent mass of H2O in H2 */
tmass_H2O += mass_H2 * (molar_mass_H2O / molar_mass_H2 );
}
if ( FLAG_CH4 ){
/* equivalent mass of H2 in CH4 */
tmass_H2O += mass_CH4 * (2 * molar_mass_H2O / molar_mass_CH4 );
}
tmass_H2O *= (*Ap->mantle_mass_ptr); /* total non-dimensional mass */
p_H2O = tmass_H2O / scaling2; /* equivalent surface pressure */
tmass_H2O *= scaling; /* total physical mass */
//tmass_H2O /= 1.4E21; /* non-dimensionalise according to ocean mass of H2O */
/* below I have taken the Olson and Sharp value for H2, and then scaled by the typical
molar mass for H2 and H2O, such that by construction, the ocean mass estimates from
both H2 and H2O are identical */
tmass_H2O /= 1.385185824553049e+21;
p_H2O *= SC->PRESSURE / 1.0E5; /* to bar */
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent present-day mass of ocean water from H2O (non-dimensional)",(double)tmass_H2O);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent atmospheric pressure of H2O (bar)",(double)p_H2O);CHKERRQ(ierr);
}
/* CO */
if( FLAG_CO ){
tmass_CO = mass_CO;
if ( FLAG_CO2 ){
/* equivalent mass of CO in CO2 */
tmass_CO += mass_CO2 * (molar_mass_CO / molar_mass_CO2 );
}
if ( FLAG_CH4 ){
/* equivalent mass of CO in CH4 */
tmass_CO += mass_CH4 * (molar_mass_CO / molar_mass_CH4 );
}
tmass_CO *= (*Ap->mantle_mass_ptr); /* total non-dimensional mass */
p_CO = tmass_CO / scaling2; /* equivalent surface pressure */
tmass_CO *= scaling; /* total physical mass */
tmass_CO /= 2.153674821914003e+21; /* non-dimensionalise according to ocean mass of CO */
p_CO *= SC->PRESSURE / 1.0E5; /* to bar */
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent present-day mass of ocean water from CO (non-dimensional)",(double)tmass_CO);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent atmospheric pressure of CO (bar)",(double)p_CO);CHKERRQ(ierr);
}
/* CO2 */
if( FLAG_CO2 ){
tmass_CO2 = mass_CO2;
if ( FLAG_CO ){
/* equivalent mass of H2O in H2 */
tmass_CO2 += mass_CO * (molar_mass_CO2 / molar_mass_CO );
}
if ( FLAG_CH4 ){
/* equivalent mass of CO in CH4 */
tmass_CO2 += mass_CH4 * (molar_mass_CO2 / molar_mass_CH4 );
}
tmass_CO2 *= (*Ap->mantle_mass_ptr); /* total non-dimensional mass */
p_CO2 = tmass_CO2 / scaling2; /* equivalent surface pressure */
tmass_CO2 *= scaling; /* total physical mass */
//tmass_H2O /= 1.4E21; /* non-dimensionalise according to ocean mass of H2O */
/* below I have taken the Olson and Sharp value for H2, and then scaled by the typical
molar mass for H2 and H2O, such that by construction, the ocean mass estimates from
both H2 and H2O are identical */
tmass_CO2 /= 3.383906780165486e+21;
p_CO2 *= SC->PRESSURE / 1.0E5; /* to bar */
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent present-day mass of ocean water from CO2 (non-dimensional)",(double)tmass_CO2);CHKERRQ(ierr);
ierr = PetscPrintf(PETSC_COMM_WORLD,"%-30s %-15.6g\n","Equivalent atmospheric pressure of CO2 (bar)",(double)p_CO2);CHKERRQ(ierr);
}
if( FLAG_H2 || FLAG_H2O || FLAG_CO || FLAG_CO2 || FLAG_CH4 ){
ierr = PetscPrintf(PETSC_COMM_WORLD,"************************************************\n\n");CHKERRQ(ierr);
}
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_atmosphere_pseudo_volatiles( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
PetscInt i;
Parameters const P = E->parameters;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
Atmosphere *A = &E->atmosphere;
PetscFunctionBeginUser;
for (i=0; i<Ap->n_volatiles; ++i) {
/* recall that here A->volatiles[i].p is log10(P(Pa)) from lookup */
SetInterp1dValue( Ap->volatile_parameters[i]->TP_interp, A->tsurf, &A->volatiles[i].p, NULL );
/* convert to non-dimensional (scaled) pressure */
A->volatiles[i].p -= PetscLog10Real( SC->PRESSURE );
A->volatiles[i].p = PetscPowScalar( 10.0, A->volatiles[i].p );
}
/* mass reaction has not been updated, so the imposed pressures must already adhere to
any imposed chemical equilibrium */
for (i=0; i<Ap->n_reactions; ++i) {
A->mass_reaction[i] = 0.0;
}
ierr = set_solution_from_partial_pressures( E, sol );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_atmosphere_from_ocean_moles( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
PetscInt i;
PetscScalar abund;
AtmosphereParameters Ap = E->parameters->atmosphere_parameters;
Parameters const P = E->parameters;
FundamentalConstants const FC = P->fundamental_constants;
ScalingConstants const SC = P->scaling_constants;
PetscFunctionBeginUser;
/* convert from ocean moles to mass and update total abundance */
for (i=0; i<Ap->n_volatiles; ++i) {
/* mass of volatile */
abund = Ap->volatile_parameters[i]->initial_ocean_moles * FC->OCEAN_MOLES;
abund *= Ap->volatile_parameters[i]->molar_mass;
abund /= (*Ap->mantle_mass_ptr) * SC->VOLATILE * 4.0 * PETSC_PI;
Ap->volatile_parameters[i]->initial_total_abundance = abund;
}
ierr = set_ic_atmosphere_from_initial_total_abundance( E, sol );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_atmosphere_from_initial_total_abundance( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
PetscFunctionBeginUser;
ierr = solve_for_initial_partial_pressure( E );CHKERRQ(ierr);
/* below will also update mass reaction terms, which can be
non-zero */
ierr = set_solution_from_partial_pressures( E, sol );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode set_ic_atmosphere_from_partial_pressure( Ctx *E, Vec sol )
{
PetscErrorCode ierr;
PetscInt i;
Atmosphere *A = &E->atmosphere;
AtmosphereParameters Ap = E->parameters->atmosphere_parameters;
PetscFunctionBeginUser;
/* set initial partial pressure to A->volatiles[i].p */
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].p = Ap->volatile_parameters[i]->initial_atmos_pressure;
}
/* mass reaction has not been updated, so the imposed pressures must already adhere to
any imposed chemical equilibrium */
for (i=0; i<Ap->n_reactions; ++i) {
A->mass_reaction[i] = 0.0;
}
ierr = set_solution_from_partial_pressures( E, sol );CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode solve_for_initial_partial_pressure( Ctx *E )
{
/* returns both volatile partial pressures and reaction
masses, both of which can be non-zero */
PetscErrorCode ierr;
SNES snes;
Vec x,r;
PetscScalar *xx;
PetscInt i;
Atmosphere *A = &E->atmosphere;
Parameters const P = E->parameters;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscFunctionBeginUser;
ierr = SNESCreate( PETSC_COMM_WORLD, &snes );CHKERRQ(ierr);
/* Use this to address this specific SNES (nonlinear solver) from the command
line or options file, e.g. -atmosic_snes_view */
ierr = SNESSetOptionsPrefix(snes,"atmosic_");CHKERRQ(ierr);
ierr = VecCreate( PETSC_COMM_WORLD, &x );CHKERRQ(ierr);
/* Vector with one dof per volatile species, and one dof per volatile reaction */
ierr = VecSetSizes( x, PETSC_DECIDE, Ap->n_volatiles + Ap->n_reactions );CHKERRQ(ierr);
ierr = VecSetFromOptions(x);CHKERRQ(ierr);
ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
ierr = SNESSetFunction(snes,r,objective_function_initial_partial_pressure,E);CHKERRQ(ierr);
/* initialise vector x with initial guess */
ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
for (i=0; i<Ap->n_volatiles; ++i) {
xx[i] = Ap->volatile_parameters[i]->initial_atmos_pressure;
}
/* initial guesses for reaction masses */
for (i=Ap->n_volatiles; i<Ap->n_volatiles + Ap->n_reactions; ++i) {
xx[i] = 0.0; /* assume we close to equilibrium */
}
ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
/* it's often tricky for the user to know how the mass is partitioned between
species, and the solver fails if the initial guess of the mass distribution
is way off. The newtontr seems to help, and particularly increasing the size
of delta0 to 10.0. But I'm sure further optimisations are possible */
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_type","newtontr");CHKERRQ(ierr);
/* Inform the nonlinear solver to generate a finite-difference approximation
to the Jacobian */
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_mf",NULL);CHKERRQ(ierr);
/* Turn off convergence based on step size */
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_stol","0");CHKERRQ(ierr);
/* Turn off convergenced based on trust region tolerance */
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_trtol","0");CHKERRQ(ierr);
/* https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESNEWTONTR.html */
/* initial size of the trust region is delta0*norm2(x) */
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_tr_delta0","10.0");CHKERRQ(ierr);
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_rtol","1.0e-9");CHKERRQ(ierr);
ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_atol","1.0e-9");CHKERRQ(ierr);
ierr = PetscOptionsSetValue(NULL,"-atmosic_ksp_rtol","1.0e-9");CHKERRQ(ierr);
ierr = PetscOptionsSetValue(NULL,"-atmosic_ksp_rtol","1.0e-9");CHKERRQ(ierr);
/* below previously suggested by PS, but not currently used */
//ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_linesearch_damping","0.01");CHKERRQ(ierr);
//ierr = PetscOptionsSetValue(NULL,"-atmosic_snes_max_it","10000");CHKERRQ(ierr);
/* For solver analysis/debugging/tuning, activate a custom monitor with a flag */
{
PetscBool flg = PETSC_FALSE;
ierr = PetscOptionsGetBool(NULL,NULL,"-atmosic_snes_verbose_monitor",&flg,NULL);CHKERRQ(ierr);
if (flg) {
ierr = SNESMonitorSet(snes,SNESMonitorVerbose,NULL,NULL);CHKERRQ(ierr);
}
}
/* Solve */
ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); /* Picks up any additional options (note prefix) */
ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
{
SNESConvergedReason reason;
ierr = SNESGetConvergedReason(snes,&reason);CHKERRQ(ierr);
if (reason < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_CONV_FAILED,
"Nonlinear solver didn't converge: %s\n",SNESConvergedReasons[reason]);
}
ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
for (i=0; i<Ap->n_volatiles; ++i) {
if( A->volatiles[i].p < 0.0 ){
/* Sanity check on solution (since it's non-unique) */
SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_CONV_FAILED,
"Unphysical initial volatile partial pressure: volatile %d, x: %g",i,A->volatiles[i].p);
}
else{
A->volatiles[i].p = xx[i];
}
}
for (i=0; i<Ap->n_reactions; ++i) {
A->mass_reaction[i] = xx[Ap->n_volatiles+i];
}
ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
ierr = VecDestroy(&x);CHKERRQ(ierr);
ierr = VecDestroy(&r);CHKERRQ(ierr);
ierr = SNESDestroy(&snes);CHKERRQ(ierr);
PetscFunctionReturn(0);
}
static PetscErrorCode objective_function_initial_partial_pressure( SNES snes, Vec x, Vec f, void *ptr)
{
PetscErrorCode ierr;
const PetscScalar *xx;
PetscScalar *ff;
PetscInt i;
Ctx *E = (Ctx*) ptr;
Atmosphere *A = &E->atmosphere;
Parameters const P = E->parameters;
FundamentalConstants const FC = P->fundamental_constants;
ScalingConstants const SC = P->scaling_constants;
AtmosphereParameters const Ap = P->atmosphere_parameters;
PetscFunctionBeginUser;
ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
ierr = VecGetArray(f,&ff);CHKERRQ(ierr);
for (i=0; i<Ap->n_volatiles; ++i) {
A->volatiles[i].p = xx[i];
}
for (i=0; i<Ap->n_reactions; ++i) {
A->mass_reaction[i] = xx[Ap->n_volatiles+i];
}
ierr = set_reservoir_volatile_content( A, Ap, FC, SC ); CHKERRQ(ierr);
/* mass conservation for each volatile */