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 "radiation.hpp"
9 #include "pellet.hpp"
10 #include "thermal_bath.hpp"
11 
12 enum class CollisionMode{
13  None = 0,
14  MonteCarloNonConserving, // monte-carlo non-conserving
15  NonLinearFokkerPlanckLandau, // fully non-linear multi-species Fokker-Planck-Landau [Hager JCP 315, 644-660 (2016)]
16  PETScLandauAMR // PETSc multi-species Landau with adaptive mesh refinement (ongoing implementation)
17 };
18 
19 struct Sources{
20  private:
21 
22  StepTrigger step_trigger; // Triggers when f_source is called
23 
24  public:
25 
27 
28  // List of sources (other collisions-related members will be folded into a child class of Source)
37 
41 
42  Sources(NLReader::NamelistReader& nlr, const DomainDecomposition<DeviceType>& pol_decomp, const MagneticField<DeviceType>& magnetic_field, const Grid<DeviceType>& grid, const VelocityGrid& vgrid, Plasma& plasma, bool overwrite_existing_files) {
43  nlr.use_namelist("sml_param");
44  f_source_period = nlr.get<int>("sml_f_source_period", 1); // Period (in time steps) of when f_source is called
45 
47 
48  // Set up for each source:
49  const int CALL_AT_F_SOURCE = 1; // Set these triggers so that these sources are called every time f_source is called
50  bool sml_radiation = nlr.get<bool>("sml_radiation",false); // Switch for a simple radiation module using part of the ADAS
51  // data base (must be configured with the rad_param namelist)
52  bool sml_neutral = nlr.get<bool>("sml_neutral",false); // Switch for using the neutrals
53  bool sml_source = nlr.get<bool>("sml_source",false); // Switch for using the heat source
54  bool sml_current_drive = nlr.get<bool>("sml_current_drive_on",false); // Switch for current drive
55  bool sml_thermal_bath_on = nlr.get<bool>("sml_thermal_bath_on",false); // Switch for thermal bath
56  if (sml_radiation) sml_neutral = true; // radiation overrides neutrals (should throw error instead)
57 
58  if(sml_radiation){
59  radiation = Radiation(nlr, grid, magnetic_field);
60  }
61 
62  // Temporary: 'f0_grid' is true if every species is either adiabatic or dynamic_f0
63  bool f0_grid = plasma.true_for_all_species([&](const Species<DeviceType>& species){return (species.is_adiabatic || species.dynamic_f0);});
64 
65  if(sml_neutral){
66  if(!f0_grid) exit_XGC("\nERROR: sml_neutrals can't (yet) be used unless all non-adiabatic species have ptl_dynamic_f0 = .true.\n");
67 
68  nlr.use_namelist("neu_param");
69  int start_time = nlr.get<int>("neu_start_time",2);
70  neutrals.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
71  }
72 
73  if(sml_source){
74  if(!f0_grid) exit_XGC("\nERROR: sml_source can't (yet) be used unless all non-adiabatic species have ptl_dynamic_f0 = .true.\n");
75 
76  nlr.use_namelist("src_param");
77  int start_time = nlr.get<int>("src_start_time",1);
78  heat_and_torque.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
79  }
80 
81  if(sml_current_drive){
82  nlr.use_namelist("src_param");
83  int start_time = nlr.get<int>("src_current_drive_start_time",1); // Time step at which current drive turns on
84  current_drive.step_trigger = StepTrigger(CALL_AT_F_SOURCE, start_time);
85  }
86 
87  if(sml_thermal_bath_on){
89  }
90 
91  nlr.use_namelist("diff_param");
92  bool diff_on = nlr.get<bool>("diff_on", false);
93  if(diff_on){
94  // If diffusion is on, call diffusion constructor
95  GPTLstart("F_DIFF_INIT");
96  diffusion = Diffusion(nlr, magnetic_field, grid);
97  GPTLstop("F_DIFF_INIT");
98  }
99 
100  nlr.use_namelist("src_param");
101  bool pellet_on = nlr.get<bool>("src_pellet_on", false);
102  if(pellet_on){
103  // If pellet is on, call pellet constructor
104  pellet = Pellet(nlr, grid, plasma, magnetic_field);
105  }
106 
107  // Collisions
108  nlr.use_namelist("col_param");
109  int col_mode_in = nlr.get<int>("col_mode", 0); // Collision method. 0: No collisions. 1: uses the monte carlo non-conserving collisions. 4: uses the non-linear Fokker Planck Landau collision operator. 5: PETSc Landau AMR collisions.
110  if(col_mode_in==0) col_mode = CollisionMode::None;
111  else if(col_mode_in==1) col_mode = CollisionMode::MonteCarloNonConserving;
112  else if(col_mode_in==4) col_mode = CollisionMode::NonLinearFokkerPlanckLandau;
113  else if(col_mode_in==5) col_mode = CollisionMode::PETScLandauAMR;
114  else{
115  if(col_mode_in==2){
116  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");
117  }
118  exit_XGC("\nInvalid col_mode specified. Acceptable options are 0, 1, 4, or 5.\n");
119  }
120 
122  col_grid = CollisionGrid<DeviceType>(nlr, magnetic_field, grid, overwrite_existing_files);
124  col_mc = MonteCarloCollider<DeviceType>(nlr, pol_decomp, magnetic_field, grid, plasma.n_nonadiabatic_species);
125  }
126 
128  collisions.step_trigger = StepTrigger(CALL_AT_F_SOURCE, 0); // start time is zero
129  }
130  }
131 
132  // Close streams for diagnostics associated with a source. This could be part of a destructor.
135  col_grid.io_stream->Close();
136  }
137  }
138 
139  /* Returns whether the sources should be executed on this step
140  * */
141  bool is_triggered(int step){
142  return step_trigger.is_triggered(step);
143  }
144 };
145 
146 #endif
Source neutrals
Definition: sources.hpp:29
Sources(NLReader::NamelistReader &nlr, const DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field, const Grid< DeviceType > &grid, const VelocityGrid &vgrid, Plasma &plasma, bool overwrite_existing_files)
Definition: sources.hpp:42
Definition: sources.hpp:19
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
bool is_rank_zero()
Definition: globals.hpp:27
bool dynamic_f0
Whether f0 can evolve in time.
Definition: species.hpp:96
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
Source collisions
Definition: sources.hpp:33
int f_source_period
Frequency of f_sources.
Definition: sources.hpp:26
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:1224
Definition: thermal_bath.hpp:8
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
Definition: pellet.hpp:12
MonteCarloCollider< DeviceType > col_mc
Monte carlo collisions, only used by MonteCarloNonConserving.
Definition: sources.hpp:39
Source current_drive
Definition: sources.hpp:34
Definition: source.hpp:7
Definition: diffusion.hpp:17
CollisionGrid< DeviceType > col_grid
Collision grid, only used by NonLinearFokkerPlanckLandau.
Definition: sources.hpp:40
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
Source heat_and_torque
Definition: sources.hpp:32
bool is_triggered(int step) const
Definition: step_trigger.hpp:25
CollisionMode
Definition: sources.hpp:12
int n_nonadiabatic_species
Number of nonadiabatic species.
Definition: plasma.hpp:90
Diffusion diffusion
Definition: sources.hpp:30
StepTrigger step_trigger
Definition: source.hpp:8
CollisionMode col_mode
Collision operator model.
Definition: sources.hpp:38
void exit_XGC(std::string msg)
Definition: globals.hpp:37
bool is_adiabatic
Whether this species is adiabatic.
Definition: species.hpp:80
Pellet pellet
Definition: sources.hpp:35
Definition: magnetic_field.F90:1
Definition: plasma.hpp:13
Radiation radiation
Definition: sources.hpp:31
ThermalBath thermal_bath
Definition: sources.hpp:36
Definition: species.hpp:75
std::shared_ptr< XGC_IO_Stream > io_stream
Definition: col_grid.hpp:214
void close_streams()
Definition: sources.hpp:133
StepTrigger step_trigger
Definition: sources.hpp:22
bool is_triggered(int step)
Definition: sources.hpp:141
Definition: step_trigger.hpp:4
bool true_for_all_species(F func) const
Definition: plasma.hpp:190
bool diag_on
Definition: col_grid.hpp:213
static int GPTLstop(const char *name)
Definition: timer_macro.hpp:10
Definition: radiation.hpp:13