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 #include "diffusion.hpp"
8 
9 enum class CollisionMode{
10  None = 0,
11  MonteCarloNonConserving, // monte-carlo non-conserving
12  NonLinearFokkerPlanckLandau, // fully non-linear multi-species Fokker-Planck-Landau [Hager JCP 315, 644-660 (2016)]
13  PETScLandauAMR // PETSc multi-species Landau with adaptive mesh refinement (ongoing implementation)
14 };
15 
16 struct Sources{
17  private:
18 
19  StepTrigger step_trigger; // Triggers when f_source is called
20 
21  public:
22 
24 
25  // List of sources (other collisions-related members will be folded into a child class of Source)
32 
36 
38  nlr.use_namelist("sml_param");
39  f_source_period = nlr.get<int>("sml_f_source_period", 1);
40 
42 
43  // Set up for each source:
44  const int CALL_AT_F_SOURCE = 1; // Set these triggers so that these sources are called every time f_source is called
45  bool sml_radiation = nlr.get<bool>("sml_radiation",false);
46  bool sml_neutral = nlr.get<bool>("sml_neutral",false);
47  bool sml_source = nlr.get<bool>("sml_source",false);
48  bool sml_current_drive = nlr.get<bool>("sml_current_drive_on",false);
49  if (sml_radiation) sml_neutral = true; // radiation overrides neutrals (should throw error instead)
50 
51  if(sml_radiation){
52  nlr.use_namelist("rad_param");
53  int start_time = nlr.get<int>("rad_start_time",1);
54  radiation.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
55  }
56 
57  if(sml_neutral){
58  nlr.use_namelist("neu_param");
59  int start_time = nlr.get<int>("neu_start_time",2);
60  neutrals.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
61  }
62 
63  if(sml_source){
64  nlr.use_namelist("src_param");
65  int start_time = nlr.get<int>("src_start_time",1);
66  heat_and_torque.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
67  }
68 
69  if(sml_current_drive){
70  nlr.use_namelist("src_param");
71  int start_time = nlr.get<int>("src_current_drive_start_time",1);
72  current_drive.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
73  }
74 
75  nlr.use_namelist("diff_param");
76  bool diff_on = nlr.get<bool>("diff_on", false);
77  if(diff_on){
78  // If diffusion is on, call diffusion constructor
79  diffusion = Diffusion(nlr, magnetic_field);
80 
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  col_grid = CollisionGrid<DeviceType>(nlr, magnetic_field);
102  col_mc = MonteCarloCollider<DeviceType>(nlr, magnetic_field, plasma.n_nonadiabatic_species);
103  }
104 
106  collisions.step_trigger = StepTrigger(CALL_AT_F_SOURCE, 0); // start time is zero
107  }
108  }
109 
110  /* Returns whether the sources should be executed on this step
111  * */
112  bool is_triggered(int step){
113  return step_trigger.is_triggered(step);
114  }
115 };
116 
117 #endif
Source neutrals
Definition: sources.hpp:26
Definition: sources.hpp:16
bool is_rank_zero()
Definition: globals.hpp:27
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
Source collisions
Definition: sources.hpp:30
int f_source_period
Frequency of f_sources.
Definition: sources.hpp:23
Definition: velocity_grid.hpp:8
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:1238
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
MonteCarloCollider< DeviceType > col_mc
Monte carlo collisions, only used by MonteCarloNonConserving.
Definition: sources.hpp:34
Sources(NLReader::NamelistReader &nlr, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid, const Plasma &plasma)
Definition: sources.hpp:37
Source current_drive
Definition: sources.hpp:31
Definition: source.hpp:7
Definition: diffusion.hpp:14
CollisionGrid< DeviceType > col_grid
Collision grid, only used by NonLinearFokkerPlanckLandau.
Definition: sources.hpp:35
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
Source heat_and_torque
Definition: sources.hpp:29
CollisionMode
Definition: sources.hpp:9
Source radiation
Definition: sources.hpp:28
int n_nonadiabatic_species
Number of nonadiabatic species.
Definition: plasma.hpp:91
bool is_triggered(int step)
Definition: step_trigger.hpp:23
Diffusion diffusion
Definition: sources.hpp:27
StepTrigger step_trigger
Definition: source.hpp:8
CollisionMode col_mode
Collision operator model.
Definition: sources.hpp:33
void exit_XGC(std::string msg)
Definition: globals.hpp:37
Definition: magnetic_field.F90:1
Definition: plasma.hpp:14
StepTrigger step_trigger
Definition: sources.hpp:19
bool is_triggered(int step)
Definition: sources.hpp:112
Definition: step_trigger.hpp:4