forked from fmihpc/vlasiator
-
Notifications
You must be signed in to change notification settings - Fork 0
/
parameters.cpp
924 lines (842 loc) · 47.2 KB
/
parameters.cpp
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
/*
* This file is part of Vlasiator.
* Copyright 2010-2016 Finnish Meteorological Institute
*
* For details of usage, see the COPYING file and read the "Rules of the Road"
* at http://www.physics.helsinki.fi/vlasiator/
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License along
* with this program; if not, write to the Free Software Foundation, Inc.,
* 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
*/
#include "parameters.h"
#include "object_wrapper.h"
#include "particle_species.h"
#include "readparameters.h"
#include <algorithm>
#include <cstdlib>
#include <iostream>
#include <limits>
#include <set>
#include <unistd.h>
#include "fieldtracing/fieldtracing.h"
#ifndef NAN
#define NAN 0
#endif
using namespace std;
typedef Parameters P;
// Using numeric_limits<Real>::max() leads to FP exceptions inside boost programoptions, use a slightly smaller value to
// avoid...
const Real LARGE_REAL = 1e20;
// Define static members:
int P::geometry = geometry::XYZ6D;
Real P::xmin = NAN;
Real P::xmax = NAN;
Real P::ymin = NAN;
Real P::ymax = NAN;
Real P::zmin = NAN;
Real P::zmax = NAN;
Real P::dx_ini = NAN;
Real P::dy_ini = NAN;
Real P::dz_ini = NAN;
uint P::xcells_ini = numeric_limits<uint>::max();
uint P::ycells_ini = numeric_limits<uint>::max();
uint P::zcells_ini = numeric_limits<uint>::max();
Real P::t = 0;
Real P::t_min = 0;
Real P::t_max = LARGE_REAL;
Real P::dt = NAN;
Real P::vlasovSolverMaxCFL = NAN;
Real P::vlasovSolverMinCFL = NAN;
Real P::fieldSolverMaxCFL = NAN;
Real P::fieldSolverMinCFL = NAN;
uint P::fieldSolverSubcycles = 1;
bool P::amrTransShortPencils = false;
uint P::tstep = 0;
uint P::tstep_min = 0;
uint P::tstep_max = 0;
uint P::diagnosticInterval = numeric_limits<uint>::max();
bool P::writeInitialState = true;
bool P::writeFullBGB = false;
bool P::meshRepartitioned = true;
bool P::prepareForRebalance = false;
vector<CellID> P::localCells;
vector<string> P::systemWriteName;
vector<string> P::systemWritePath;
vector<Real> P::systemWriteTimeInterval;
vector<int> P::systemWriteDistributionWriteStride;
vector<int> P::systemWriteDistributionWriteXlineStride;
vector<int> P::systemWriteDistributionWriteYlineStride;
vector<int> P::systemWriteDistributionWriteZlineStride;
vector<Real> P::systemWriteDistributionWriteShellRadius;
vector<int> P::systemWriteDistributionWriteShellStride;
vector<bool> P::systemWriteFsGrid;
bool P::systemWriteAllDROs;
bool P::diagnosticWriteAllDROs;
vector<int> P::systemWrites;
vector<pair<string, string>> P::systemWriteHints;
vector<pair<string, string>> P::restartWriteHints;
vector<pair<string, string>> P::restartReadHints;
Real P::saveRestartWalltimeInterval = -1.0;
uint P::exitAfterRestarts = numeric_limits<uint>::max();
uint64_t P::vlsvBufferSize = 0;
int P::restartStripeFactor = 0;
int P::systemStripeFactor = 0;
string P::restartWritePath = string("");
uint P::transmit = 0;
bool P::recalculateStencils = true;
bool P::propagateVlasovAcceleration = true;
bool P::propagateVlasovTranslation = true;
bool P::propagateField = true;
bool P::dynamicTimestep = true;
Real P::maxWaveVelocity = 0.0;
uint P::maxFieldSolverSubcycles = 0.0;
int P::maxSlAccelerationSubcycles = 0.0;
Real P::resistivity = NAN;
bool P::fieldSolverDiffusiveEterms = true;
uint P::ohmHallTerm = 0;
uint P::ohmGradPeTerm = 0;
Real P::electronTemperature = 0.0;
Real P::electronDensity = 0.0;
Real P::electronPTindex = 1.0;
string P::restartFileName = string("");
bool P::isRestart = false;
int P::writeAsFloat = false;
int P::writeRestartAsFloat = false;
string P::loadBalanceAlgorithm = string("");
std::map<std::string, std::string> P::loadBalanceOptions;
uint P::rebalanceInterval = numeric_limits<uint>::max();
vector<string> P::outputVariableList;
vector<string> P::diagnosticVariableList;
string P::projectName = string("");
bool P::vlasovAccelerateMaxwellianBoundaries = false;
Real P::maxSlAccelerationRotation = 10.0;
Real P::hallMinimumRhom = physicalconstants::MASS_PROTON;
Real P::hallMinimumRhoq = physicalconstants::CHARGE;
bool P::bailout_write_restart = false;
Real P::bailout_min_dt = NAN;
Real P::bailout_max_memory = 1073741824.;
uint P::bailout_velocity_space_wall_margin = 0;
uint P::vamrMaxVelocityRefLevel = 0;
Realf P::vamrRefineLimit = 1.0;
Realf P::vamrCoarsenLimit = 0.5;
string P::vamrVelRefCriterion = string("");
int P::amrMaxSpatialRefLevel = 0;
int P::amrMaxAllowedSpatialRefLevel = -1;
bool P::adaptRefinement = false;
bool P::refineOnRestart = false;
bool P::forceRefinement = false;
bool P::shouldFilter = false;
bool P::useAlpha = true;
Real P::alphaRefineThreshold = 0.5;
Real P::alphaCoarsenThreshold = -1.0;
bool P::useJPerB = true;
Real P::jperbRefineThreshold = 0.5;
Real P::jperbCoarsenThreshold = -1.0;
Real P::alphaDRhoWeight = 1.0;
Real P::alphaDUWeight = 1.0;
Real P::alphaDPSqWeight = 1.0;
Real P::alphaDBSqWeight = 1.0;
Real P::alphaDBWeight = 1.0;
uint P::refineCadence = 5;
Real P::refineAfter = 0.0;
Real P::refineRadius = LARGE_REAL;
int P::maxFilteringPasses = 0;
uint P::amrBoxHalfWidthX = 1;
uint P::amrBoxHalfWidthY = 1;
uint P::amrBoxHalfWidthZ = 1;
Realf P::amrBoxCenterX = 0.0;
Realf P::amrBoxCenterY = 0.0;
Realf P::amrBoxCenterZ = 0.0;
vector<string> P::blurPassString;
std::vector<int> P::numPasses; //numpasses
std::string tracerString; /*!< Fieldline tracer to use for coupling ionosphere and magnetosphere */
bool P::computeCurvature;
bool P::addParameters() {
typedef Readparameters RP;
// the other default parameters we read through the add/get interface
RP::add("io.diagnostic_write_interval", "Write diagnostic output every arg time steps", numeric_limits<uint>::max());
RP::addComposing(
"io.system_write_t_interval",
"Save the simulation every arg simulated seconds. Negative values disable writes. [Define for all groups.]");
RP::addComposing("io.system_write_file_name",
"Save the simulation to this file name series. [Define for all groups.]");
RP::addComposing("io.system_write_path",
"Save this series in this location. Default is ./ [Define for all groups or none.]");
RP::addComposing("io.system_write_distribution_stride",
"Every this many cells write out their velocity space. 0 is none. [Define for all groups.]");
RP::addComposing("io.system_write_distribution_xline_stride",
"Every this many lines of cells along the x direction write out their velocity space. 0 is none. "
"[Define for all groups.]");
RP::addComposing("io.system_write_distribution_yline_stride",
"Every this many lines of cells along the y direction write out their velocity space. 0 is none. "
"[Define for all groups.]");
RP::addComposing("io.system_write_distribution_zline_stride",
"Every this many lines of cells along the z direction write out their velocity space. 0 is none. "
"[Define for all groups.]");
RP::addComposing("io.system_write_distribution_shell_radius",
"At cells intersecting spheres with those radii centred at the origin write out their velocity "
"space. 0 is none.");
RP::addComposing("io.system_write_distribution_shell_stride",
"Every this many cells for those on selected shells write out their velocity space. 0 is none.");
RP::addComposing("io.system_write_fsgrid_variables", "If 0 don't write fsgrid DROs, if 1 do write them.");
RP::addComposing(
"io.system_write_mpiio_hint_key",
"MPI-IO hint key passed to the non-restart IO. Has to be matched by io.system_write_mpiio_hint_value.");
RP::addComposing(
"io.system_write_mpiio_hint_value",
"MPI-IO hint value passed to the non-restart IO. Has to be matched by io.system_write_mpiio_hint_key.");
RP::addComposing(
"io.restart_write_mpiio_hint_key",
"MPI-IO hint key passed to the restart IO. Has to be matched by io.restart_write_mpiio_hint_value.");
RP::addComposing(
"io.restart_write_mpiio_hint_value",
"MPI-IO hint value passed to the restart IO. Has to be matched by io.restart_write_mpiio_hint_key.");
RP::addComposing(
"io.restart_read_mpiio_hint_key",
"MPI-IO hint key passed to the restart IO. Has to be matched by io.restart_read_mpiio_hint_value.");
RP::addComposing(
"io.restart_read_mpiio_hint_value",
"MPI-IO hint value passed to the restart IO. Has to be matched by io.restart_read_mpiio_hint_key.");
RP::add("io.write_initial_state",
"Write initial state, not even the 0.5 dt propagation is done. Do not use for restarting. ", false);
RP::add("io.write_full_bgb_data", "Write a dedicated file containing all BGB components and first derivatives, then exit.", false);
RP::add("io.restart_walltime_interval",
"Save the complete simulation in given walltime intervals. Negative values disable writes.", -1.0);
RP::add("io.number_of_restarts", "Exit the simulation after certain number of walltime-based restarts.",
numeric_limits<uint>::max());
RP::add("io.vlsv_buffer_size",
"Buffer size passed to VLSV writer (bytes, up to uint64_t), default 0 as this is sensible on sisu", 0);
RP::add("io.write_restart_stripe_factor", "Stripe factor for restart and initial grid writing. Default 0 to inherit.", 0);
RP::add("io.write_system_stripe_factor", "Stripe factor for bulk file writing. Default 0 to inherit.", 0);
RP::add("io.write_as_float", "If true, write in floats instead of doubles", false);
RP::add("io.restart_write_path",
"Path to the location where restart files should be written. Defaults to the local directory, also if the "
"specified destination is not writeable.",
string("./"));
RP::add("propagate_field", "Propagate magnetic field during the simulation", true);
RP::add("propagate_vlasov_acceleration",
"Propagate distribution functions during the simulation in velocity space. If false, it is propagated with "
"zero length timesteps.",
true);
RP::add("propagate_vlasov_translation",
"Propagate distribution functions during the simulation in ordinary space. If false, it is propagated with "
"zero length timesteps.",
true);
RP::add("dynamic_timestep", "If true, timestep is set based on CFL limits (default on)", true);
RP::add("hallMinimumRho",
"Minimum rho value used for the Hall and electron pressure gradient terms in the Lorentz force and in the "
"field solver. Default is very low and has no effect in practice.",
1.0);
RP::add("project",
"Specify the name of the project to use. Supported to date (20150610): Alfven Diffusion Dispersion "
"Distributions Firehose Flowthrough Fluctuations Harris KHB Larmor Magnetosphere Multipeak Riemann1 Shock "
"Shocktest Template test_fp testHall test_trans VelocityBox verificationLarmor",
string(""));
RP::add("restart.write_as_float", "If true, write restart fields in floats instead of doubles", false);
RP::add("restart.filename", "Restart from this vlsv file. No restart if empty file.", string(""));
RP::add("gridbuilder.geometry", "Simulation geometry XY4D,XZ4D,XY5D,XZ5D,XYZ6D", string("XYZ6D"));
RP::add("gridbuilder.x_min", "Minimum value of the x-coordinate.", NAN);
RP::add("gridbuilder.x_max", "Minimum value of the x-coordinate.", NAN);
RP::add("gridbuilder.y_min", "Minimum value of the y-coordinate.", NAN);
RP::add("gridbuilder.y_max", "Minimum value of the y-coordinate.", NAN);
RP::add("gridbuilder.z_min", "Minimum value of the z-coordinate.", NAN);
RP::add("gridbuilder.z_max", "Minimum value of the z-coordinate.", NAN);
RP::add("gridbuilder.x_length", "Number of cells in x-direction in initial grid.", 0);
RP::add("gridbuilder.y_length", "Number of cells in y-direction in initial grid.", 0);
RP::add("gridbuilder.z_length", "Number of cells in z-direction in initial grid.", 0);
RP::add("gridbuilder.dt", "Initial timestep in seconds.", 0.0);
RP::add("gridbuilder.t_max",
"Maximum simulation time, in seconds. If timestep_max limit is hit first this time will never be reached",
LARGE_REAL);
RP::add("gridbuilder.timestep_max",
"Max. value for timesteps. If t_max limit is hit first, this step will never be reached",
numeric_limits<uint>::max());
// Field solver parameters
RP::add("fieldsolver.maxWaveVelocity",
"Maximum wave velocity allowed in the fastest velocity determination in m/s, default unlimited", LARGE_REAL);
RP::add("fieldsolver.maxSubcycles", "Maximum allowed field solver subcycles", 1);
RP::add("fieldsolver.resistivity", "Resistivity for the eta*J term in Ohm's law.", 0.0);
RP::add("fieldsolver.diffusiveEterms", "Enable diffusive terms in the computation of E", true);
RP::add(
"fieldsolver.ohmHallTerm",
"Enable/choose spatial order of the Hall term in Ohm's law. 0: off, 1: 1st spatial order, 2: 2nd spatial order",
0);
RP::add(
"fieldsolver.ohmGradPeTerm",
"Enable/choose spatial order of the electron pressure gradient term in Ohm's law. 0: off, 1: 1st spatial order.",
0);
RP::add("fieldsolver.electronTemperature",
"Upstream electron temperature to be used for the electron pressure gradient term (K).", 0.0);
RP::add("fieldsolver.electronDensity",
"Upstream electron density to be used for the electron pressure gradient term (m^-3).", 0.0);
RP::add("fieldsolver.electronPTindex",
"Polytropic index for electron pressure gradient term. 0 is isobaric, 1 is isothermal, 1.667 is adiabatic "
"electrons, ",
0.0);
RP::add("fieldsolver.maxCFL",
"The maximum CFL limit for field propagation. Used to set timestep if dynamic_timestep is true.", 0.5);
RP::add("fieldsolver.minCFL",
"The minimum CFL limit for field propagation. Used to set timestep if dynamic_timestep is true.", 0.4);
// Vlasov solver parameters
RP::add("vlasovsolver.maxSlAccelerationRotation",
"Maximum rotation angle (degrees) allowed by the Semi-Lagrangian solver (Use >25 values with care)", 25.0);
RP::add("vlasovsolver.maxSlAccelerationSubcycles", "Maximum number of subcycles for acceleration", 1);
RP::add("vlasovsolver.maxCFL",
"The maximum CFL limit for vlasov propagation in ordinary space. Used to set timestep if dynamic_timestep "
"is true.",
0.99);
RP::add("vlasovsolver.minCFL",
"The minimum CFL limit for vlasov propagation in ordinary space. Used to set timestep if dynamic_timestep "
"is true.",
0.8);
RP::add("vlasovsolver.accelerateMaxwellianBoundaries",
"Propagate maxwellian boundary cell contents in velocity space. Default false.",
false);
// Load balancing parameters
RP::add("loadBalance.algorithm", "Load balancing algorithm to be used", string("RCB"));
RP::add("loadBalance.tolerance", "Load imbalance tolerance", string("1.05"));
RP::add("loadBalance.rebalanceInterval", "Load rebalance interval (steps)", 10);
RP::addComposing("loadBalance.optionKey", "Zoltan option key. Has to be matched by loadBalance.optionValue.");
RP::addComposing("loadBalance.optionValue", "Zoltan option value. Has to be matched by loadBalance.optionKey.");
// Output variable parameters
RP::add("io.system_write_all_data_reducers", "If 0 don't write all DROs, if 1 do write them.", false);
// NOTE Do not remove the : before the list of variable names as this is parsed by tools/check_vlasiator_cfg.sh
RP::addComposing("variables.output",
string() +
"List of data reduction operators (DROs) to add to the grid file output. Each variable to be "
"added has to be on a new line output = XXX. Names are case insensitive. " +
"Available (20230628): " + "fg_b fg_b_background fg_b_perturbed fg_b_background_vol fg_derivs_b_background fg_e " +
"vg_rhom vg_rhoq populations_vg_rho " + "fg_rhom fg_rhoq " + "vg_v fg_v populations_vg_v " +
"populations_vg_moments_thermal populations_vg_moments_nonthermal " +
"populations_vg_effectivesparsitythreshold populations_vg_rho_loss_adjust " +
"populations_vg_energydensity populations_vg_precipitationdifferentialflux " +
"populations_vg_heatflux " +
"populations_vg_nonmaxwellianity " +
"vg_maxdt_acceleration vg_maxdt_translation populations_vg_maxdt_acceleration " +
"populations_vg_maxdt_translation " +
"fg_maxdt_fieldsolver " + "vg_rank fg_rank fg_amr_level vg_loadbalance_weight " +
"vg_boundarytype fg_boundarytype vg_boundarylayer fg_boundarylayer " +
"populations_vg_blocks vg_f_saved " + "populations_vg_acceleration_subcycles " +
"vg_e_vol fg_e_vol " +
"fg_e_hall vg_e_gradpe fg_b_vol vg_b_vol vg_b_background_vol vg_b_perturbed_vol " +
"vg_pressure fg_pressure populations_vg_ptensor " + "vg_b_vol_derivatives fg_derivs " +
"ig_fac ig_latitude ig_chi0 ig_cellarea ig_upmappedarea ig_sigmap ig_sigmah ig_sigmaparallel ig_rhon " +
"ig_electrontemp ig_solverinternals ig_upmappednodecoords ig_upmappedb ig_openclosed ig_potential "+
"ig_precipitation ig_deltaphi "+
"ig_inplanecurrent ig_b ig_e vg_drift vg_ionospherecoupling vg_connection vg_fluxrope fg_curvature "+
"vg_amr_drho vg_amr_du vg_amr_dpsq vg_amr_dbsq vg_amr_db vg_amr_alpha vg_amr_reflevel vg_amr_jperb "+
"vg_amr_translate_comm vg_gridcoordinates fg_gridcoordinates ");
RP::addComposing(
"variables_deprecated.output",
string() + "List of deprecated names for data reduction operators (DROs). Names are case insensitive. " +
"Available (20190521): " + "B BackgroundB fg_BackgroundB PerturbedB fg_PerturbedB " + "E " +
"Rhom Rhoq populations_Rho " + "V populations_V " +
"populations_moments_Backstream populations_moments_NonBackstream " +
"populations_moments_thermal populations_moments_nonthermal " +
"populations_minvalue populations_EffectiveSparsityThreshold populations_RhoLossAdjust "
"populations_rho_loss_adjust" +
"populations_EnergyDensity populations_PrecipitationFlux populations_precipitationdifferentialflux" +
"LBweight vg_lbweight vg_loadbalanceweight MaxVdt MaxRdt populations_MaxVdt populations_MaxRdt " +
"populations_maxdt_acceleration populations_maxdt_translation MaxFieldsdt fg_maxfieldsdt" +
"MPIrank FsGridRank " + "FsGridBoundaryType BoundaryType FsGridBoundaryLayer BoundaryLayer " +
"populations_Blocks fSaved vg_fsaved" + "populations_accSubcycles populations_acceleration_subcycles" +
"VolE vg_VolE Evol E_vol fg_VolE fg_Evol " +
"HallE fg_HallE GradPeE e_gradpe VolB vg_VolB fg_VolB B_vol Bvol vg_Bvol fg_volB fg_Bvol " +
"BackgroundVolB PerturbedVolB " + "Pressure vg_Pressure fg_Pressure populations_PTensor " +
"BVOLderivs b_vol_derivs");
RP::add("io.diagnostic_write_all_data_reducers", "Write all available diagnostic reducers", false);
// NOTE Do not remove the : before the list of variable names as this is parsed by tools/check_vlasiator_cfg.sh
RP::addComposing("variables.diagnostic",
string() +
"List of data reduction operators (DROs) to add to the diagnostic runtime output. Each "
"variable to be added has to be on a new line diagnostic = XXX. Names are case insensitive. " +
"Available (20221221): " + "populations_vg_blocks " +
"vg_rhom populations_vg_rho_loss_adjust " + "vg_loadbalance_weight " +
"vg_maxdt_acceleration vg_maxdt_translation " + "fg_maxdt_fieldsolver " +
"populations_vg_maxdt_acceleration populations_vg_maxdt_translation " +
"populations_vg_maxdistributionfunction populations_vg_mindistributionfunction");
RP::addComposing("variables_deprecated.diagnostic",
string() +
"List of deprecated data reduction operators (DROs) to add to the diagnostic runtime output. "
"Names are case insensitive. " +
"Available (20201111): " + "rhom populations_rholossadjust populations_rho_loss_adjust " +
"populations_blocks lbweight loadbalance_weight " + "vg_lbweight vg_loadbalanceweight " +
"maxvdt maxdt_acceleration " + "maxrdt maxdt_translation " +
"populations_maxvdt populations_maxrdt " +
"populations_maxdt_acceleration populations_maxdt_translation " +
"populations_maxdistributionfunction populations_mindistributionfunction " +
"maxfieldsdt maxdt_fieldsolver fg_maxfieldsdt");
// bailout parameters
RP::add("bailout.write_restart",
"If 1, write a restart file on bailout. Gets reset when sending a STOP (1) or a KILL (0).", true);
RP::add("bailout.min_dt", "Minimum time step below which bailout occurs (s).", 1e-6);
RP::add("bailout.max_memory", "Maximum amount of memory used per node (in GiB) over which bailout occurs.",
1073741824.);
RP::add("bailout.velocity_space_wall_block_margin", "Distance from the velocity space limits in blocks, if the distribution function reaches that distance from the wall we bail out to avoid hitting the wall.", 1);
// Velocity Refinement parameters
RP::add("VAMR.vel_refinement_criterion", "Name of the velocity refinement criterion", string(""));
RP::add("VAMR.max_velocity_level", "Maximum velocity mesh refinement level", (uint)0);
RP::add("VAMR.refine_limit",
"If the refinement criterion function returns a larger value than this, block is refined", (Realf)1.0);
RP::add("VAMR.coarsen_limit",
"If the refinement criterion function returns a smaller value than this, block can be coarsened",
(Realf)0.5);
// Spatial Refinement parameters
RP::add("AMR.max_spatial_level", "Maximum absolute spatial mesh refinement level", (uint)0);
RP::add("AMR.max_allowed_spatial_level", "Maximum currently allowed spatial mesh refinement level", -1);
RP::add("AMR.should_refine","If false, do not refine Vlasov grid regardless of max spatial level",true);
RP::add("AMR.adapt_refinement","If true, re-refine vlasov grid every refine_cadence balance", false);
RP::add("AMR.refine_on_restart","If true, re-refine vlasov grid on restart. DEPRECATED, consider using the DOMR command", false);
RP::add("AMR.force_refinement","If true, refine/unrefine the vlasov grid to match the config on restart", false);
RP::add("AMR.should_filter","If true, filter vlasov grid with boxcar filter on restart",false);
RP::add("AMR.use_alpha1","Use the maximum of dimensionless gradients alpha_1 as a refinement index", true);
RP::add("AMR.alpha1_refine_threshold","Determines the minimum value of alpha_1 to refine cells", 0.5);
RP::add("AMR.alpha1_coarsen_threshold","Determines the maximum value of alpha_1 to unrefine cells, default half of the refine threshold", -1.0);
RP::add("AMR.use_alpha2","Use J/B_perp as a refinement index", true);
RP::add("AMR.alpha2_refine_threshold","Determines the minimum value of alpha_2 to refine cells", 0.5);
RP::add("AMR.alpha2_coarsen_threshold","Determines the maximum value of alpha_2 to unrefine cells, default half of the refine threshold", -1.0);
RP::add("AMR.refine_cadence","Refine every nth load balance", 5);
RP::add("AMR.refine_after","Start refinement after this many simulation seconds", 0.0);
RP::add("AMR.refine_radius","Maximum distance from Earth to refine", LARGE_REAL);
RP::add("AMR.alpha1_drho_weight","Multiplier for delta rho in alpha calculation", 1.0);
RP::add("AMR.alpha1_du_weight","Multiplier for delta U in alpha calculation", 1.0);
RP::add("AMR.alpha1_dpsq_weight","Multiplier for delta p squared in alpha calculation", 1.0);
RP::add("AMR.alpha1_dbsq_weight","Multiplier for delta B squared in alpha calculation", 1.0);
RP::add("AMR.alpha1_db_weight","Multiplier for delta B in alpha calculation", 1.0);
RP::add("AMR.box_half_width_x", "Half width of the box that is refined (for testing)", (uint)1);
RP::add("AMR.box_half_width_y", "Half width of the box that is refined (for testing)", (uint)1);
RP::add("AMR.box_half_width_z", "Half width of the box that is refined (for testing)", (uint)1);
RP::add("AMR.box_center_x", "x coordinate of the center of the box that is refined (for testing)", 0.0);
RP::add("AMR.box_center_y", "y coordinate of the center of the box that is refined (for testing)", 0.0);
RP::add("AMR.box_center_z", "z coordinate of the center of the box that is refined (for testing)", 0.0);
RP::add("AMR.transShortPencils", "if true, use one-cell pencils", false);
RP::addComposing("AMR.filterpasses", string("AMR filter passes for each individual refinement level"));
RP::add("fieldtracing.fieldLineTracer", "Field line tracing method to use for coupling ionosphere and magnetosphere (options are: Euler, BS)", std::string("Euler"));
RP::add("fieldtracing.tracer_max_allowed_error", "Maximum allowed error for the adaptive field line tracers ", 1000);
RP::add("fieldtracing.tracer_max_attempts", "Maximum allowed attempts for the adaptive field line tracers", 100);
RP::add("fieldtracing.tracer_min_dx", "Minimum allowed field line tracer step length for the adaptive field line tracers (m)", 100e3);
RP::add("fieldtracing.fullbox_and_fluxrope_max_absolute_distance_to_trace", "Maximum absolute distance in m to trace along the field line before ending. Defaults to the sum of the simulation box edge lengths LX+LY+LZ if set <= 0.", -1);
RP::add("fieldtracing.fullbox_max_incomplete_cells", "Maximum fraction of cells left incomplete when stopping tracing loop for full box tracing. Defaults to zero to process all, will be slow at scale! Both fluxrope_max_incomplete_cells and fullbox_max_incomplete_cells will be achieved.", 0);
RP::add("fieldtracing.fluxrope_max_incomplete_cells", "Maximum fraction of cells left incomplete when stopping loop for flux rope tracing. Defaults to zero to process all, will be slow at scale! Both fluxrope_max_incomplete_cells and fullbox_max_incomplete_cells will be achieved.", 0);
RP::add("fieldtracing.use_reconstruction_cache", "Use the cache to store reconstruction coefficients. (0: don't, 1: use)", 0);
RP::add("fieldtracing.fluxrope_max_curvature_radii_to_trace", "Maximum number of seedpoint curvature radii to trace forward and backward from each DCCRG cell to find flux ropes", 10);
RP::add("fieldtracing.fluxrope_max_curvature_radii_extent", "Maximum extent in seedpoint curvature radii from the seed a field line is allowed to extend to be counted as a flux rope", 2);
return true;
}
void Parameters::getParameters() {
typedef Readparameters RP;
// get numerical values of the parameters
RP::get("io.diagnostic_write_interval", P::diagnosticInterval);
RP::get("io.diagnostic_write_all_data_reducers", P::diagnosticWriteAllDROs);
RP::get("io.system_write_t_interval", P::systemWriteTimeInterval);
RP::get("io.system_write_file_name", P::systemWriteName);
RP::get("io.system_write_path", P::systemWritePath);
RP::get("io.system_write_distribution_stride", P::systemWriteDistributionWriteStride);
RP::get("io.system_write_distribution_xline_stride", P::systemWriteDistributionWriteXlineStride);
RP::get("io.system_write_distribution_yline_stride", P::systemWriteDistributionWriteYlineStride);
RP::get("io.system_write_distribution_zline_stride", P::systemWriteDistributionWriteZlineStride);
RP::get("io.system_write_distribution_shell_radius", P::systemWriteDistributionWriteShellRadius);
RP::get("io.system_write_distribution_shell_stride", P::systemWriteDistributionWriteShellStride);
RP::get("io.system_write_fsgrid_variables", P::systemWriteFsGrid);
RP::get("io.system_write_all_data_reducers", P::systemWriteAllDROs);
RP::get("io.write_initial_state", P::writeInitialState);
RP::get("io.write_full_bgb_data", P::writeFullBGB);
RP::get("io.restart_walltime_interval", P::saveRestartWalltimeInterval);
RP::get("io.number_of_restarts", P::exitAfterRestarts);
RP::get("io.vlsv_buffer_size", P::vlsvBufferSize);
RP::get("io.write_restart_stripe_factor", P::restartStripeFactor);
RP::get("io.write_system_stripe_factor", P::systemStripeFactor);
RP::get("io.restart_write_path", P::restartWritePath);
RP::get("io.write_as_float", P::writeAsFloat);
// Checks for validity of io and restart parameters
int myRank;
MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
const string prefix = string("./");
if (access(&(P::restartWritePath[0]), W_OK) != 0) {
if (myRank == MASTER_RANK) {
cerr << "ERROR restart write path " << P::restartWritePath << " not writeable, defaulting to local directory."
<< endl;
}
P::restartWritePath = prefix;
}
size_t maxSize = 0;
maxSize = max(maxSize, P::systemWriteTimeInterval.size());
maxSize = max(maxSize, P::systemWriteName.size());
maxSize = max(maxSize, P::systemWritePath.size());
maxSize = max(maxSize, P::systemWriteDistributionWriteStride.size());
maxSize = max(maxSize, P::systemWriteDistributionWriteXlineStride.size());
maxSize = max(maxSize, P::systemWriteDistributionWriteYlineStride.size());
maxSize = max(maxSize, P::systemWriteDistributionWriteZlineStride.size());
if (P::systemWriteTimeInterval.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_t_interval should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteName.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_file_name should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWritePath.size() != maxSize && P::systemWritePath.size() != 0) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_path should be defined for all file types or none at all." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteDistributionWriteStride.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_distribution_stride should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteDistributionWriteXlineStride.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_distribution_xline_stride should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteDistributionWriteYlineStride.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_distribution_yline_stride should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteDistributionWriteZlineStride.size() != maxSize) {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_distribution_zline_stride should be defined for all file types." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteDistributionWriteShellStride.size() != P::systemWriteDistributionWriteShellRadius.size()) {
if (myRank == MASTER_RANK) {
cerr << "ERROR You should set the same number of io.system_write_distribution_shell_stride "
<< "and io.system_write_distribution_shell_radius." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
if (P::systemWriteFsGrid.size() != maxSize) {
if (P::systemWriteFsGrid.size() == 0) {
for (uint i = 0; i < maxSize; i++) {
P::systemWriteFsGrid.push_back(true);
}
} else {
if (myRank == MASTER_RANK) {
cerr << "ERROR io.system_write_fsgrid_variables should be defined for all file types (or none at all)." << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
}
}
if (P::systemWritePath.size() == 0) {
for (uint i = 0; i < P::systemWriteName.size(); i++) {
P::systemWritePath.push_back(string("./"));
}
} else {
for (uint i = 0; i < P::systemWritePath.size(); i++) {
if (access(&(P::systemWritePath.at(i)[0]), W_OK) != 0) {
if (myRank == MASTER_RANK) {
cerr << "ERROR " << P::systemWriteName.at(i) << " write path " << P::systemWritePath.at(i)
<< " not writeable, defaulting to local directory." << endl;
}
P::systemWritePath.at(i) = prefix;
}
}
}
vector<string> mpiioKeys, mpiioValues;
RP::get("io.system_write_mpiio_hint_key", mpiioKeys);
RP::get("io.system_write_mpiio_hint_value", mpiioValues);
if (mpiioKeys.size() != mpiioValues.size()) {
if (myRank == MASTER_RANK) {
cerr << "WARNING the number of io.system_write_mpiio_hint_key and io.system_write_mpiio_hint_value do not "
"match. Disregarding these options."
<< endl;
}
} else {
for (uint i = 0; i < mpiioKeys.size(); i++) {
P::systemWriteHints.push_back({mpiioKeys[i], mpiioValues[i]});
}
}
mpiioKeys.clear();
mpiioValues.clear();
RP::get("io.restart_write_mpiio_hint_key", mpiioKeys);
RP::get("io.restart_write_mpiio_hint_value", mpiioValues);
if (mpiioKeys.size() != mpiioValues.size()) {
if (myRank == MASTER_RANK) {
cerr << "WARNING the number of io.restart_write_mpiio_hint_key and io.restart_write_mpiio_hint_value do not "
"match. Disregarding these options."
<< endl;
}
} else {
for (uint i = 0; i < mpiioKeys.size(); i++) {
P::restartWriteHints.push_back({mpiioKeys[i], mpiioValues[i]});
}
}
mpiioKeys.clear();
mpiioValues.clear();
RP::get("io.restart_read_mpiio_hint_key", mpiioKeys);
RP::get("io.restart_read_mpiio_hint_value", mpiioValues);
if (mpiioKeys.size() != mpiioValues.size()) {
if (myRank == MASTER_RANK) {
cerr << "WARNING the number of io.restart_read_mpiio_hint_key and io.restart_read_mpiio_hint_value do not "
"match. Disregarding these options."
<< endl;
}
} else {
for (uint i = 0; i < mpiioKeys.size(); i++) {
P::restartReadHints.push_back({mpiioKeys[i], mpiioValues[i]});
}
}
RP::get("propagate_field", P::propagateField);
RP::get("propagate_vlasov_acceleration", P::propagateVlasovAcceleration);
RP::get("propagate_vlasov_translation", P::propagateVlasovTranslation);
RP::get("dynamic_timestep", P::dynamicTimestep);
Real hallRho;
RP::get("hallMinimumRho", hallRho);
P::hallMinimumRhom = hallRho * physicalconstants::MASS_PROTON;
P::hallMinimumRhoq = hallRho * physicalconstants::CHARGE;
RP::get("restart.write_as_float", P::writeRestartAsFloat);
RP::get("restart.filename", P::restartFileName);
P::isRestart = (P::restartFileName != string(""));
RP::get("project", P::projectName);
if (RP::helpRequested) {
P::projectName = string("Magnetosphere");
}
/*get numerical values, let Readparameters handle the conversions*/
string geometryString;
RP::get("gridbuilder.geometry", geometryString);
RP::get("gridbuilder.x_min", P::xmin);
RP::get("gridbuilder.x_max", P::xmax);
RP::get("gridbuilder.y_min", P::ymin);
RP::get("gridbuilder.y_max", P::ymax);
RP::get("gridbuilder.z_min", P::zmin);
RP::get("gridbuilder.z_max", P::zmax);
RP::get("gridbuilder.x_length", P::xcells_ini);
RP::get("gridbuilder.y_length", P::ycells_ini);
RP::get("gridbuilder.z_length", P::zcells_ini);
RP::get("VAMR.max_velocity_level", P::vamrMaxVelocityRefLevel);
RP::get("VAMR.vel_refinement_criterion", P::vamrVelRefCriterion);
RP::get("VAMR.refine_limit", P::vamrRefineLimit);
RP::get("VAMR.coarsen_limit", P::vamrCoarsenLimit);
RP::get("AMR.max_spatial_level", P::amrMaxSpatialRefLevel);
RP::get("AMR.max_allowed_spatial_level", P::amrMaxAllowedSpatialRefLevel);
if(P::amrMaxAllowedSpatialRefLevel < 0) { // negative (default is -1) just goes to max
P::amrMaxAllowedSpatialRefLevel = P::amrMaxSpatialRefLevel; // set max allowed to the same as the absolute max
}
if(P::amrMaxSpatialRefLevel < P::amrMaxAllowedSpatialRefLevel) {
if(myRank == MASTER_RANK) {
cerr << "AMR.max_allowed_spatial_level cannot be greater than AMR.max_spatial_level!\n";
}
MPI_Abort(MPI_COMM_WORLD, 1);
}
RP::get("AMR.adapt_refinement",P::adaptRefinement);
RP::get("AMR.refine_on_restart",P::refineOnRestart);
RP::get("AMR.force_refinement",P::forceRefinement);
RP::get("AMR.should_filter",P::shouldFilter);
RP::get("AMR.use_alpha1",P::useAlpha);
RP::get("AMR.alpha1_refine_threshold",P::alphaRefineThreshold);
RP::get("AMR.alpha1_coarsen_threshold",P::alphaCoarsenThreshold);
if (P::useAlpha && P::alphaCoarsenThreshold < 0) {
P::alphaCoarsenThreshold = P::alphaRefineThreshold / 2.0;
}
if (P::useAlpha && P::alphaRefineThreshold < 0) {
if (myRank == MASTER_RANK) {
cerr << "ERROR invalid alpha_1 refine threshold" << endl;
}
MPI_Abort(MPI_COMM_WORLD, 1);
}
RP::get("AMR.use_alpha2",P::useJPerB);
RP::get("AMR.alpha2_refine_threshold",P::jperbRefineThreshold);
RP::get("AMR.alpha2_coarsen_threshold",P::jperbCoarsenThreshold);
if (P::useJPerB && P::jperbCoarsenThreshold < 0) {
P::jperbCoarsenThreshold = P::jperbRefineThreshold / 2.0;
}
if (P::useJPerB && P::jperbRefineThreshold < 0) {
if (myRank == MASTER_RANK) {
cerr << "ERROR invalid alpha_2 refine threshold" << endl;
}
MPI_Abort(MPI_COMM_WORLD, 1);
}
RP::get("AMR.refine_cadence",P::refineCadence);
RP::get("AMR.refine_after",P::refineAfter);
RP::get("AMR.refine_radius",P::refineRadius);
RP::get("AMR.alpha1_drho_weight", P::alphaDRhoWeight);
RP::get("AMR.alpha1_du_weight", P::alphaDUWeight);
RP::get("AMR.alpha1_dpsq_weight", P::alphaDPSqWeight);
RP::get("AMR.alpha1_dbsq_weight", P::alphaDBSqWeight);
RP::get("AMR.alpha1_db_weight", P::alphaDBWeight);
RP::get("AMR.box_half_width_x", P::amrBoxHalfWidthX);
RP::get("AMR.box_half_width_y", P::amrBoxHalfWidthY);
RP::get("AMR.box_half_width_z", P::amrBoxHalfWidthZ);
RP::get("AMR.box_center_x", P::amrBoxCenterX);
RP::get("AMR.box_center_y", P::amrBoxCenterY);
RP::get("AMR.box_center_z", P::amrBoxCenterZ);
RP::get("AMR.transShortPencils", P::amrTransShortPencils);
RP::get("AMR.filterpasses", P::blurPassString);
// If we are in an AMR run we need to set up the filtering scheme.
if (P::amrMaxSpatialRefLevel>0){
bool isEmpty = blurPassString.size() == 0;
if (!isEmpty){
//sanity check=> user should define a pass for every level
if ((int)blurPassString.size() != P::amrMaxSpatialRefLevel + 1) {
cerr << "Filter Passes=" << blurPassString.size() << "\t" << "AMR Levels=" << P::amrMaxSpatialRefLevel + 1 << endl;
cerr << "FilterPasses do not match AMR levels. \t" << " in " << __FILE__ << ":" << __LINE__ << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
//sort the filtering passes per refLevel
numPasses.clear();
//Parse to a vector of ints
for (auto pass : blurPassString){
P::numPasses.push_back(stoi(pass));
}
sort(numPasses.begin(),numPasses.end(),greater<int>());
}else{
//here we will default to manually constructing the number of passes
numPasses.clear();
auto g_sequence=[](int size){
int retval=1;
while(size!=0){
retval*=2;
size-=1;
}
return retval;
};
int maxPasses=g_sequence(P::amrMaxSpatialRefLevel-1);
for (int refLevel=0; refLevel<=P::amrMaxSpatialRefLevel; refLevel++){
numPasses.push_back(maxPasses);
maxPasses/=2;
}
//Overwrite passes for the highest refLevel. We do not want to filter there.
numPasses.at(P::amrMaxSpatialRefLevel) = 0;
}
P::maxFilteringPasses = numPasses[0];
}
if (geometryString == "XY4D") {
P::geometry = geometry::XY4D;
} else if (geometryString == "XZ4D") {
P::geometry = geometry::XZ4D;
} else if (geometryString == "XY5D") {
P::geometry = geometry::XY5D;
} else if (geometryString == "XZ5D") {
P::geometry = geometry::XZ5D;
} else if (geometryString == "XYZ6D") {
P::geometry = geometry::XYZ6D;
} else {
cerr << "Unknown simulation geometry " << geometryString << " in " << __FILE__ << ":" << __LINE__ << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
if (P::vamrCoarsenLimit >= P::vamrRefineLimit) {
cerr << "vamrRefineLimit must be smaller than vamrCoarsenLimit!" << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
if (P::xmax < P::xmin || (P::ymax < P::ymin || P::zmax < P::zmin)) {
cerr << "Box domain error!" << endl;
MPI_Abort(MPI_COMM_WORLD, 1);
}
// Set some parameter values.
P::dx_ini = (P::xmax - P::xmin) / P::xcells_ini;
P::dy_ini = (P::ymax - P::ymin) / P::ycells_ini;
P::dz_ini = (P::zmax - P::zmin) / P::zcells_ini;
RP::get("gridbuilder.dt", P::dt);
RP::get("gridbuilder.t_max", P::t_max);
RP::get("gridbuilder.timestep_max", P::tstep_max);
if (P::dynamicTimestep)
P::dt = 0.0; // if dynamic timestep then first dt is always 0
// if we are restarting, t,t_min, tstep, tstep_min will be overwritten in readGrid
P::t_min = 0;
P::t = P::t_min;
P::tstep_min = 0;
P::tstep = P::tstep_min;
// Get field solver parameters
RP::get("fieldsolver.maxWaveVelocity", P::maxWaveVelocity);
RP::get("fieldsolver.maxSubcycles", P::maxFieldSolverSubcycles);
RP::get("fieldsolver.resistivity", P::resistivity);
RP::get("fieldsolver.diffusiveEterms", P::fieldSolverDiffusiveEterms);
RP::get("fieldsolver.ohmHallTerm", P::ohmHallTerm);
RP::get("fieldsolver.ohmGradPeTerm", P::ohmGradPeTerm);
RP::get("fieldsolver.electronTemperature", P::electronTemperature);
RP::get("fieldsolver.electronDensity", P::electronDensity);
RP::get("fieldsolver.electronPTindex", P::electronPTindex);
RP::get("fieldsolver.maxCFL", P::fieldSolverMaxCFL);
RP::get("fieldsolver.minCFL", P::fieldSolverMinCFL);
// Get Vlasov solver parameters
RP::get("vlasovsolver.maxSlAccelerationRotation", P::maxSlAccelerationRotation);
RP::get("vlasovsolver.maxSlAccelerationSubcycles", P::maxSlAccelerationSubcycles);
RP::get("vlasovsolver.maxCFL", P::vlasovSolverMaxCFL);
RP::get("vlasovsolver.minCFL", P::vlasovSolverMinCFL);
RP::get("vlasovsolver.accelerateMaxwellianBoundaries", P::vlasovAccelerateMaxwellianBoundaries);
// Get load balance parameters
RP::get("loadBalance.algorithm", P::loadBalanceAlgorithm);
loadBalanceOptions["IMBALANCE_TOL"] = "";
RP::get("loadBalance.tolerance", loadBalanceOptions["IMBALANCE_TOL"]);
RP::get("loadBalance.rebalanceInterval", P::rebalanceInterval);
std::vector<std::string> loadBalanceKeys;
std::vector<std::string> loadBalanceValues;
RP::get("loadBalance.optionKey", loadBalanceKeys);
RP::get("loadBalance.optionValue", loadBalanceValues);
if (loadBalanceKeys.size() != loadBalanceValues.size()) {
if (myRank == MASTER_RANK) {
cerr << "WARNING the number of load balance keys and values do not match. Disregarding these options." << endl;
}
} else {
for (size_t i = 0; i < loadBalanceKeys.size(); ++i) {
loadBalanceOptions[loadBalanceKeys[i]] = loadBalanceValues[i];
}
}
// Get output variable parameters
RP::get("variables.output", P::outputVariableList);
RP::get("variables.diagnostic", P::diagnosticVariableList);
// Filter duplicate variable names
set<string> dummy(P::outputVariableList.begin(), P::outputVariableList.end());
P::outputVariableList.clear();
P::outputVariableList.insert(P::outputVariableList.end(), dummy.begin(), dummy.end());
dummy.clear();
dummy.insert(P::diagnosticVariableList.begin(), P::diagnosticVariableList.end());
P::diagnosticVariableList.clear();
P::diagnosticVariableList.insert(P::diagnosticVariableList.end(), dummy.begin(), dummy.end());
// Get parameters related to bailout
RP::get("bailout.write_restart", P::bailout_write_restart);
RP::get("bailout.min_dt", P::bailout_min_dt);
RP::get("bailout.max_memory", P::bailout_max_memory);
RP::get("bailout.velocity_space_wall_block_margin", P::bailout_velocity_space_wall_margin);
if(P::bailout_velocity_space_wall_margin > MAX_BLOCKS_PER_DIM / 2 && myRank == MASTER_RANK) {
std::cerr << "bailout.velocity_space_wall_block_margin is larger than 0.5 * MAX_BLOCKS_PER_DIM, aborting." << std::endl;
abort();
}
for (size_t s = 0; s < P::systemWriteName.size(); ++s) {
P::systemWrites.push_back(0);
}
RP::get("fieldtracing.fieldLineTracer", tracerString);
RP::get("fieldtracing.tracer_max_allowed_error", FieldTracing::fieldTracingParameters.max_allowed_error);
RP::get("fieldtracing.tracer_max_attempts", FieldTracing::fieldTracingParameters.max_field_tracer_attempts);
RP::get("fieldtracing.tracer_min_dx", FieldTracing::fieldTracingParameters.min_tracer_dx);
RP::get("fieldtracing.fullbox_max_incomplete_cells", FieldTracing::fieldTracingParameters.fullbox_max_incomplete_cells);
RP::get("fieldtracing.fluxrope_max_incomplete_cells", FieldTracing::fieldTracingParameters.fluxrope_max_incomplete_cells);
RP::get("fieldtracing.fullbox_and_fluxrope_max_absolute_distance_to_trace", FieldTracing::fieldTracingParameters.fullbox_and_fluxrope_max_distance);
RP::get("fieldtracing.use_reconstruction_cache", FieldTracing::fieldTracingParameters.useCache);
RP::get("fieldtracing.fluxrope_max_curvature_radii_to_trace", FieldTracing::fieldTracingParameters.fluxrope_max_curvature_radii_to_trace);
RP::get("fieldtracing.fluxrope_max_curvature_radii_extent", FieldTracing::fieldTracingParameters.fluxrope_max_curvature_radii_extent);
if(tracerString == "Euler") {
FieldTracing::fieldTracingParameters.tracingMethod = FieldTracing::Euler;
} else if (tracerString == "ADPT_Euler") {
FieldTracing::fieldTracingParameters.tracingMethod = FieldTracing::ADPT_Euler;
} else if (tracerString == "BS") {
FieldTracing::fieldTracingParameters.tracingMethod = FieldTracing::BS;
} else if (tracerString == "DP") {
FieldTracing::fieldTracingParameters.tracingMethod = FieldTracing::DPrince;
} else {
cerr << __FILE__ << ":" << __LINE__ << " ERROR: Unknown value for fieldtracing.fieldLineTracer: " << tracerString << endl;
abort();
}
}