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 #include "thermal_bath.hpp"
9 
10 enum class CollisionMode{
11  None = 0,
12  MonteCarloNonConserving, // monte-carlo non-conserving
13  NonLinearFokkerPlanckLandau, // fully non-linear multi-species Fokker-Planck-Landau [Hager JCP 315, 644-660 (2016)]
14  PETScLandauAMR // PETSc multi-species Landau with adaptive mesh refinement (ongoing implementation)
15 };
16 
17 struct Sources{
18  private:
19 
20  StepTrigger step_trigger; // Triggers when f_source is called
21 
22  public:
23 
25 
26  // List of sources (other collisions-related members will be folded into a child class of Source)
34 
38 
39  Sources(NLReader::NamelistReader& nlr, const MagneticField<DeviceType>& magnetic_field, const Grid<DeviceType>& grid, const VelocityGrid& vgrid, const Plasma& plasma, bool overwrite_existing_files) {
40  nlr.use_namelist("sml_param");
41  f_source_period = nlr.get<int>("sml_f_source_period", 1);
42 
44 
45  // Set up for each source:
46  const int CALL_AT_F_SOURCE = 1; // Set these triggers so that these sources are called every time f_source is called
47  bool sml_radiation = nlr.get<bool>("sml_radiation",false);
48  bool sml_neutral = nlr.get<bool>("sml_neutral",false);
49  bool sml_source = nlr.get<bool>("sml_source",false);
50  bool sml_current_drive = nlr.get<bool>("sml_current_drive_on",false);
51  bool sml_thermal_bath_on = nlr.get<bool>("sml_thermal_bath_on",false);
52  if (sml_radiation) sml_neutral = true; // radiation overrides neutrals (should throw error instead)
53 
54  if(sml_radiation){
55  nlr.use_namelist("rad_param");
56  int start_time = nlr.get<int>("rad_start_time",1);
57  radiation.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
58  }
59 
60  if(sml_neutral){
61  nlr.use_namelist("neu_param");
62  int start_time = nlr.get<int>("neu_start_time",2);
63  neutrals.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
64  }
65 
66  if(sml_source){
67  nlr.use_namelist("src_param");
68  int start_time = nlr.get<int>("src_start_time",1);
69  heat_and_torque.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
70  }
71 
72  if(sml_current_drive){
73  nlr.use_namelist("src_param");
74  int start_time = nlr.get<int>("src_current_drive_start_time",1);
75  current_drive.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
76  }
77 
78  if(sml_thermal_bath_on){
80  }
81 
82  nlr.use_namelist("diff_param");
83  bool diff_on = nlr.get<bool>("diff_on", false);
84  if(diff_on){
85  // If diffusion is on, call diffusion constructor
86  diffusion = Diffusion(nlr, magnetic_field);
87 
88  int start_time = nlr.get<int>("diff_start_time", 1);
89  diffusion.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
90  }
91 
92  // Collisions
93  nlr.use_namelist("col_param");
94  int col_mode_in = nlr.get<int>("col_mode", 0);
95  if(col_mode_in==0) col_mode = CollisionMode::None;
96  else if(col_mode_in==1) col_mode = CollisionMode::MonteCarloNonConserving;
97  else if(col_mode_in==4) col_mode = CollisionMode::NonLinearFokkerPlanckLandau;
98  else if(col_mode_in==5) col_mode = CollisionMode::PETScLandauAMR;
99  else{
100  if(col_mode_in==2){
101  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");
102  }
103  exit_XGC("\nInvalid col_mode specified. Acceptable options are 0, 1, 4, or 5.\n");
104  }
105 
107  col_grid = CollisionGrid<DeviceType>(nlr, magnetic_field, overwrite_existing_files);
109  col_mc = MonteCarloCollider<DeviceType>(nlr, magnetic_field, grid, plasma.n_nonadiabatic_species);
110  }
111 
113  collisions.step_trigger = StepTrigger(CALL_AT_F_SOURCE, 0); // start time is zero
114  }
115  }
116 
117  // Close streams for diagnostics associated with a source. This could be part of a destructor.
120  col_grid.io_stream->Close();
121  }
122  }
123 
124  /* Returns whether the sources should be executed on this step
125  * */
126  bool is_triggered(int step){
127  return step_trigger.is_triggered(step);
128  }
129 };
130 
131 #endif
Source neutrals
Definition: sources.hpp:27
Definition: sources.hpp:17
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:31
int f_source_period
Frequency of f_sources.
Definition: sources.hpp:24
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:1230
Definition: thermal_bath.hpp:8
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
Sources(NLReader::NamelistReader &nlr, const MagneticField< DeviceType > &magnetic_field, const Grid< DeviceType > &grid, const VelocityGrid &vgrid, const Plasma &plasma, bool overwrite_existing_files)
Definition: sources.hpp:39
MonteCarloCollider< DeviceType > col_mc
Monte carlo collisions, only used by MonteCarloNonConserving.
Definition: sources.hpp:36
Source current_drive
Definition: sources.hpp:32
Definition: source.hpp:7
Definition: diffusion.hpp:15
CollisionGrid< DeviceType > col_grid
Collision grid, only used by NonLinearFokkerPlanckLandau.
Definition: sources.hpp:37
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
Source heat_and_torque
Definition: sources.hpp:30
bool is_triggered(int step) const
Definition: step_trigger.hpp:23
CollisionMode
Definition: sources.hpp:10
Source radiation
Definition: sources.hpp:29
int n_nonadiabatic_species
Number of nonadiabatic species.
Definition: plasma.hpp:92
Diffusion diffusion
Definition: sources.hpp:28
StepTrigger step_trigger
Definition: source.hpp:8
CollisionMode col_mode
Collision operator model.
Definition: sources.hpp:35
void exit_XGC(std::string msg)
Definition: globals.hpp:37
Definition: magnetic_field.F90:1
Definition: plasma.hpp:14
ThermalBath thermal_bath
Definition: sources.hpp:33
std::shared_ptr< XGC_IO_Stream > io_stream
Definition: col_grid.hpp:211
void close_streams()
Definition: sources.hpp:118
StepTrigger step_trigger
Definition: sources.hpp:20
bool is_triggered(int step)
Definition: sources.hpp:126
Definition: step_trigger.hpp:4
bool diag_on
Definition: col_grid.hpp:210