XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
sources.hpp
Go to the documentation of this file.
1 #ifndef SOURCES_HPP
2 #define SOURCES_HPP
3 
4 #include "step_trigger.hpp"
6 #include "col_grid.hpp"
7 
8 enum class CollisionMode{
9  None = 0,
10  MonteCarloNonConserving, // monte-carlo non-conserving
11  NonLinearFokkerPlanckLandau, // fully non-linear multi-species Fokker-Planck-Landau [Hager JCP 315, 644-660 (2016)]
12  PETScLandauAMR // PETSc multi-species Landau with adaptive mesh refinement (ongoing implementation)
13 };
14 
15 // This will be a parent class for all sources
16 // For now, just contains the step triggers
17 struct Source{
18  StepTrigger step_trigger; // Triggers when it is time to call this source
19 
20  /* Returns whether this source should be executed on this step
21  * */
22  bool is_triggered(int step){
23  return step_trigger.is_triggered(step);
24  }
25 };
26 
27 struct Sources{
28  private:
29 
30  StepTrigger step_trigger; // Triggers when f_source is called
31 
32  public:
33 
35 
36  // List of sources (other collisions-related members will be folded into a child class of Source)
42 
46 
48  nlr.use_namelist("sml_param");
49  f_source_period = nlr.get<int>("sml_f_source_period", 1);
50 
52 
53  // Set up for each source:
54  const int CALL_AT_F_SOURCE = 1; // Set these triggers so that these sources are called every time f_source is called
55  bool sml_radiation = nlr.get<bool>("sml_radiation",false);
56  bool sml_neutral = nlr.get<bool>("sml_neutral",false);
57  bool sml_source = nlr.get<bool>("sml_source",false);
58  if (sml_radiation) sml_neutral = true; // radiation overrides neutrals (should throw error instead)
59 
60  if(sml_radiation){
61  nlr.use_namelist("rad_param");
62  int start_time = nlr.get<int>("rad_start_time",1);
63  radiation.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
64  }
65 
66  if(sml_neutral){
67  nlr.use_namelist("neu_param");
68  int start_time = nlr.get<int>("neu_start_time",2);
69  neutrals.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
70  }
71 
72  if(sml_source){
73  nlr.use_namelist("src_param");
74  int start_time = nlr.get<int>("src_start_time",1);
75  heat_and_torque.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
76  }
77 
78  nlr.use_namelist("diff_param");
79  bool diff_on = nlr.get<bool>("diff_on", false);
80  if(diff_on){
81  int start_time = nlr.get<int>("diff_start_time", 1);
82  diffusion.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
83  }
84 
85  // Collisions
86  nlr.use_namelist("col_param");
87  int col_mode_in = nlr.get<int>("col_mode", 0);
88  if(col_mode_in==0) col_mode = CollisionMode::None;
89  else if(col_mode_in==1) col_mode = CollisionMode::MonteCarloNonConserving;
90  else if(col_mode_in==4) col_mode = CollisionMode::NonLinearFokkerPlanckLandau;
91  else if(col_mode_in==5) col_mode = CollisionMode::PETScLandauAMR;
92  else{
93  if(col_mode_in==2){
94  if(is_rank_zero()) printf("\ncol_mode==2 is currently not supported. If you are reviving it, dont forget to address this potential issue: https://github.com/PrincetonUniversity/XGC-Devel/issues/266");
95  }
96  exit_XGC("\nInvalid col_mode specified. Acceptable options are 0, 1, 4, or 5.\n");
97  }
98 
100 #ifndef USE_FORTRAN_COLLISIONS
102 #endif
104  col_mc = MonteCarloCollider<DeviceType>(nlr, magnetic_field, plasma.n_nonadiabatic_species);
105  }
106 
108  collisions.step_trigger = StepTrigger(CALL_AT_F_SOURCE, 0); // start time is zero
109  }
110  }
111 
112  /* Returns whether the sources should be executed on this step
113  * */
114  bool is_triggered(int step){
115  return step_trigger.is_triggered(step);
116  }
117 };
118 
119 #endif
Source neutrals
Definition: sources.hpp:37
Definition: sources.hpp:27
bool is_rank_zero()
Definition: globals.hpp:26
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:372
Source collisions
Definition: sources.hpp:41
int f_source_period
Frequency of f_sources.
Definition: sources.hpp:34
Definition: velocity_grid.hpp:7
subroutine plasma(grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out)
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using pl...
Definition: neutral_totalf.F90:1235
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
bool is_triggered(int step)
Definition: sources.hpp:22
MonteCarloCollider< DeviceType > col_mc
Monte carlo collisions, only used by MonteCarloNonConserving.
Definition: sources.hpp:44
Sources(NLReader::NamelistReader &nlr, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, const Plasma &plasma)
Definition: sources.hpp:47
Definition: sources.hpp:17
CollisionGrid< DeviceType > col_grid
Collision grid, only used by NonLinearFokkerPlanckLandau.
Definition: sources.hpp:45
Source diffusion
Definition: sources.hpp:38
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:354
Source heat_and_torque
Definition: sources.hpp:40
CollisionMode
Definition: sources.hpp:8
Source radiation
Definition: sources.hpp:39
int n_nonadiabatic_species
Number of nonadiabatic species.
Definition: plasma.hpp:84
bool is_triggered(int step)
Definition: step_trigger.hpp:21
StepTrigger step_trigger
Definition: sources.hpp:18
CollisionMode col_mode
Collision operator model.
Definition: sources.hpp:43
void exit_XGC(std::string msg)
Definition: globals.hpp:36
Definition: magnetic_field.F90:1
Definition: plasma.hpp:15
StepTrigger step_trigger
Definition: sources.hpp:30
bool is_triggered(int step)
Definition: sources.hpp:114
Definition: step_trigger.hpp:4