XGC1
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 #include "current_drive.hpp"
13 
14 enum class CollisionMode{
15  None = 0,
16  MonteCarloNonConserving, // monte-carlo non-conserving
17  NonLinearFokkerPlanckLandau, // fully non-linear multi-species Fokker-Planck-Landau [Hager JCP 315, 644-660 (2016)]
18  PETScLandauAMR // PETSc multi-species Landau with adaptive mesh refinement (ongoing implementation)
19 };
20 
21 struct Sources{
22  private:
23 
24  StepTrigger step_trigger; // Triggers when f_source is called
25 
26  public:
27 
29 
30  // List of sources (other collisions-related members will be folded into a child class of Source)
40 
44 
45  Sources(){}
46 
47  Sources(NLReader::NamelistReader& nlr, const DomainDecomposition<DeviceType>& pol_decomp, const MagneticField<DeviceType>& magnetic_field, const Grid<DeviceType>& grid, const VelocityGrid& vgrid, Plasma& plasma, bool exclude_private_region, bool overwrite_existing_files);
48 
49  // Close streams for diagnostics associated with a source. This could be part of a destructor.
50  void close_streams(){
52  col_grid.io_stream->Close();
53  }
54  }
55 
56  /* Returns whether the sources should be executed on this step
57  * */
58  bool is_triggered(int step){
59  return step_trigger.is_triggered(step);
60  }
61 };
62 
63 #endif
bool diag_on
Definition: col_grid.hpp:139
std::shared_ptr< XGC_IO_Stream > io_stream
Definition: col_grid.hpp:140
Definition: current_drive.hpp:9
Implements an anomalous transport (advection-diffusion) model for the kinetic electrons in the plasma...
Definition: diffusion.hpp:43
Definition: fgrid_coarse_graining.hpp:39
Definition: magnetic_field.hpp:12
Definition: NamelistReader.hpp:193
Definition: pellet.hpp:13
Definition: plasma.hpp:13
Definition: radiation.hpp:14
Definition: step_trigger.hpp:4
bool is_triggered(int step) const
Definition: step_trigger.hpp:25
Move a fraction of the background grid distribution function (f0g) to the marker particles.
Definition: magnetic_field.F90:1
subroutine plasma(grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out, ignore_vacuum)
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using pl...
Definition: neutral_totalf.F90:1548
CollisionMode
Definition: sources.hpp:14
@ NonLinearFokkerPlanckLandau
Definition: source.hpp:7
Definition: sources.hpp:21
CollisionGrid< DeviceType > col_grid
Collision grid, only used by NonLinearFokkerPlanckLandau.
Definition: sources.hpp:43
Source collisions
Definition: sources.hpp:35
Radiation radiation
Definition: sources.hpp:33
Source heat_and_torque
Definition: sources.hpp:34
FGridCoarseGraining fgrid_coarse_graining
Definition: sources.hpp:39
Pellet pellet
Definition: sources.hpp:37
Sources()
Definition: sources.hpp:45
CollisionMode col_mode
Collision operator model.
Definition: sources.hpp:41
Source neutrals
Definition: sources.hpp:31
MonteCarloCollider< DeviceType > col_mc
Monte carlo collisions, only used by MonteCarloNonConserving.
Definition: sources.hpp:42
bool is_triggered(int step)
Definition: sources.hpp:58
CurrentDrive current_drive
Definition: sources.hpp:36
int f_source_period
Frequency of f_sources.
Definition: sources.hpp:28
StepTrigger step_trigger
Definition: sources.hpp:24
Diffusion diffusion
Definition: sources.hpp:32
ThermalBath thermal_bath
Definition: sources.hpp:38
void close_streams()
Definition: sources.hpp:50
Definition: thermal_bath.hpp:8
Definition: velocity_grid.hpp:8