chombo-discharge
Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
McPhoto Class Reference

Radiative tranfer equation solver using Monte-Carlo simulation. More...

#include <CD_McPhoto.H>

Inheritance diagram for McPhoto:
Inheritance graph
[legend]
Collaboration diagram for McPhoto:
Collaboration graph
[legend]

Public Types

enum class  WhichContainer {
  Photons , Bulk , EB , Domain ,
  Source
}
 Enum class for identifying various containers. Only used for interface reasons.
 

Public Member Functions

 McPhoto ()
 Constructor.
 
 McPhoto (const McPhoto &a_other)
 Disallowed copy constructor.
 
 McPhoto (const McPhoto &&a_other)
 Disallowed move constructor.
 
McPhotooperator= (const McPhoto &a_other)
 Disallowed assignment operator.
 
McPhotooperator= (const McPhoto &&a_other)
 Disallowed move assignement operator.
 
virtual ~McPhoto ()
 Destructor.
 
virtual bool advance (const Real a_dt, EBAMRCellData &a_phi, const EBAMRCellData &a_source, const bool a_zeroPhi=false) override
 Advance RTE and deposit photon particles on the mesh. More...
 
virtual bool isInstantaneous ()
 Instantaneous solver or not.
 
virtual void advancePhotonsInstantaneous (ParticleContainer< Photon > &a_bulkPhotons, ParticleContainer< Photon > &a_ebPhotons, ParticleContainer< Photon > &a_domainPhotons, ParticleContainer< Photon > &a_photons)
 Move photons and absorb them on various objects. More...
 
virtual void advancePhotonsTransient (ParticleContainer< Photon > &a_bulkPhotons, ParticleContainer< Photon > &a_ebPhotons, ParticleContainer< Photon > &a_domainPhotons, ParticleContainer< Photon > &a_photons, const Real a_dt)
 Move photons and absorb them on various objects. More...
 
virtual void computeNumPhysicalPhotons (EBAMRCellData &a_numPhysPhotonsTotal, EBAMRCellData &a_numPhysPhotonsPerPacket, const EBAMRCellData &a_source, const Real a_dt) const noexcept
 Compute the number of physical photons in each grid cell. More...
 
virtual void generateComputationalPhotons (ParticleContainer< Photon > &a_photons, const EBAMRCellData &a_numPhysicalPhotons, const size_t a_maxPhotonsPerCell) const noexcept
 Generate computational photons. More...
 
virtual void dirtySamplePhotons (ParticleContainer< PointParticle > &a_photons, EBAMRCellData &a_phi, const EBAMRCellData &a_numPhysicalPhotons, const size_t a_maxPhotonsPerCell) const noexcept
 Dirty-sampling method for photons. More...
 
virtual void remap ()
 Remap computational particles. This remaps m_photons.
 
virtual void remap (ParticleContainer< Photon > &a_photons)
 Remap computational particles. More...
 
virtual void depositPhotons ()
 Deposit photons on the mesh. More...
 
template<class P , const Real &(P::*)() const particleScalarField>
void depositPhotons (EBAMRCellData &a_phi, ParticleContainer< P > &a_particles, const DepositionType &a_deposition) const noexcept
 Deposit photons. More...
 
virtual void sortPhotonsByCell (const WhichContainer &a_which)
 Sort container by cell. More...
 
virtual void sortPhotonsByPatch (const WhichContainer &a_which)
 Sort container by patch. More...
 
virtual void parseOptions () override
 Parse class options.
 
virtual void parseRuntimeOptions () override
 Parse runtime options.
 
virtual void allocate () override
 Allocate internal storage.
 
virtual void preRegrid (const int a_base, const int a_oldFinestLevel) override
 preRegrid operations More...
 
virtual void deallocate () override
 Deallocate internal storage.
 
virtual void regrid (const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) override
 Regrid function for this class. More...
 
virtual void registerOperators () override
 Register operators that this solver needs.
 
virtual void computeBoundaryFlux (EBAMRIVData &a_ebFlux, const EBAMRCellData &a_phi) override
 Compute the boundary flux. More...
 
virtual void computeDomainFlux (EBAMRIFData &a_domainFlux, const EBAMRCellData &a_phi) override
 Compute the domain flux. More...
 
virtual void computeFlux (EBAMRCellData &a_flux, const EBAMRCellData &a_phi) override
 Compute the flux. More...
 
virtual void computeDensity (EBAMRCellData &a_isotropic, const EBAMRCellData &a_phi) override
 Compute isotropic radiative density from mesh solution. More...
 
virtual void clear (const WhichContainer &a_which)
 Clear data holder. More...
 
virtual void clear ()
 Clear data holder.
 
virtual void clear (ParticleContainer< Photon > &a_photons)
 Clear data holder. More...
 
virtual void clear (AMRParticles< Photon > &a_photons)
 Clear data holder. More...
 
virtual void writePlotFile () override
 Write plot file. More...
 
virtual void writePlotData (LevelData< EBCellFAB > &a_output, int &a_icomp, const std::string a_outputRealm, const int a_level) const noexcept override
 Write plot data. More...
 
virtual Vector< std::string > getPlotVariableNames () const override
 Write checkpoint data into handle. More...
 
virtual int getNumberOfPlotVariables () const override
 Get number of output variables. More...
 
virtual int countPhotons (const AMRParticles< Photon > &a_photons) const
 Count number of photons in particle list. More...
 
virtual ParticleContainer< Photon > & getPhotons ()
 Get m_photons. More...
 
virtual ParticleContainer< Photon > & getBulkPhotons ()
 Get bulk photons, i.e. photons absorbed on the mesh. More...
 
virtual ParticleContainer< Photon > & getEbPhotons ()
 Get eb Photons, i.e. photons absorbed on the EB. More...
 
virtual ParticleContainer< Photon > & getDomainPhotons ()
 Get domain photons, i.e. photons absorbed on domain edges/faces. More...
 
virtual ParticleContainer< Photon > & getSourcePhotons ()
 Get source photons. More...
 
virtual int getMaxPhotonsPerCell () const noexcept
 Get maximum number of photons generated per cell.
 
virtual int getNumSamplingPackets () const noexcept
 Get maximum number of photons generated per cell.
 
virtual void computeLoads (Vector< long long > &a_loads, const DisjointBoxLayout &a_dbl, const int a_level) const noexcept override
 Get computational loads for a specific grid level. This computes the number of photons in the bulk data holder. More...
 
- Public Member Functions inherited from RtSolver
 RtSolver ()
 Constructor.
 
virtual ~RtSolver ()
 Constructor (does nothing)
 
virtual std::string getName ()
 Get solver name.
 
virtual const std::string getRealm () const
 Get the realm where the solver lives.
 
virtual bool advance (const Real a_dt, const bool a_zeroPhi=false)
 Advance equation one time step. More...
 
virtual bool advance (const Real a_dt, EBAMRCellData &a_phi, const bool a_zeroPhi=false)
 Advance method. Advances one time step. More...
 
virtual void setRealm (const std::string a_realm)
 Set realm where this solver lives. More...
 
virtual void setRtSpecies (const RefCountedPtr< RtSpecies > &a_species)
 Set the radiative transfer species (RtSpecies) More...
 
virtual void setComputationalGeometry (const RefCountedPtr< ComputationalGeometry > a_computationalGeometry)
 Set computational geometry. More...
 
virtual void setAmr (const RefCountedPtr< AmrMesh > &a_amr)
 Set the amr object. More...
 
virtual void setPhase (phase::which_phase a_phase=phase::gas)
 Set phase. More...
 
virtual void setVerbosity (const int a_verbosity)
 Set verbosity. More...
 
virtual void setTime (const int a_step, const Real a_time, const Real a_dt)
 Set the time for this solver. More...
 
virtual void setStationary (const bool a_stationary)
 Set stationary solver or not. More...
 
virtual void sanityCheck ()
 Sanity check.
 
virtual bool isStationary ()
 Check if solver is stationary.
 
virtual void initialData ()
 Fill solver with initial data. By default, this sets internal data to zero. More...
 
virtual void setSource (const EBAMRCellData &a_source)
 Set source term. More...
 
virtual void setSource (const Real a_source)
 Set source. More...
 
virtual void setSource (const std::function< Real(const RealVect a_pos)> a_source)
 Set source. More...
 
virtual Real getTime () const
 Get current time. More...
 
virtual phase::which_phase getPhase ()
 Get the RTE phase. More...
 
virtual EBAMRCellDatagetPhi ()
 Get solver state. More...
 
virtual EBAMRCellDatagetSource ()
 Get multifluid source. More...
 
virtual EBAMRFluxDatagetKappa ()
 Get the absorption length. More...
 
virtual EBAMRIVDatagetKappaEb ()
 Get the absorption coefficient on irregular EB faces. More...
 
virtual RefCountedPtr< RtSpecies > & getSpecies ()
 Get species.
 

Protected Types

enum class  PhotonGeneration { Deterministic , Stochastic }
 Enum for interpreting how photons are generated when using fluid codes.
 
enum class  SourceType { Number , PerVol , PerVolSecond , PerSecond }
 Enum for adding flexibility in what the fluid source term contains.
 
enum class  IntersectionEB { Raycast , Bisection }
 An enum for switching between various types of EB intersection algorithms when intersecting photons with the EB. More...
 

Protected Member Functions

size_t drawPhotons (const Real a_source, const Real a_volume, const Real a_dt) const noexcept
 Draw photons in a cell and volume. More...
 
int domainBcMap (const int a_dir, const Side::LoHiSide a_side)
 Mapping function for domain boundary conditions. More...
 
Real randomExponential (const Real a_mean) const noexcept
 Random exponential trial. More...
 
template<class P , const Real &(P::*)() const particleScalarField>
void depositKappaConservative (EBAMRCellData &a_phi, ParticleContainer< P > &a_particles, const DepositionType a_deposition, const CoarseFineDeposition a_coarseFineDeposition) const noexcept
 This computes the "conservative" deposition, multiplied by kappa. More...
 
void depositNonConservative (EBAMRIVData &a_depositionNC, const EBAMRCellData &a_depositionKappaC) const noexcept
 Make the "non-conservative" kappa deposition. More...
 
void depositHybrid (EBAMRCellData &a_depositionH, EBAMRIVData &a_massDifference, const EBAMRIVData &a_depositionNC) const noexcept
 Make the hybrid deposition. Also compute the mass difference. More...
 
void parseTransparentBoundaries ()
 Turn on/off transparent boundaries.
 
void parseDirtySampling ()
 Parse dirty sampling for photons. More...
 
void parseDivergenceComputation ()
 Parse the divergence computation, i.e. if we blend with non-conservative divergence or not.
 
void parsePseudoPhotons ()
 Parse pseudophotons, i.e. parse how many photons can be generated per cell.
 
void parsePhotoGeneration ()
 Parse photogeneration type.
 
void parseInstantaneous ()
 Parse whether instantaneous propagation or not.
 
void parseSourceType ()
 Parse source term type.
 
void parseIntersectionEB ()
 Parse EB intersection algorithm.
 
void parseDeposition ()
 Parse deposition method.
 
void parsePlotVariables ()
 Parse plot variables.
 
virtual void depositPhotonsNGP (LevelData< EBCellFAB > &a_output, const ParticleContainer< Photon > &a_particles, const int a_level) const noexcept
 Do an NGP deposit on a specific grid level. Used for IO. More...
 
- Protected Member Functions inherited from RtSolver
void setEbIndexSpace (const RefCountedPtr< EBIndexSpace > &a_ebis)
 Set ebis.
 
void parseVerbosity () noexcept
 Parse verbosity.
 
virtual void writeData (LevelData< EBCellFAB > &a_output, int &a_comp, const EBAMRCellData &a_data, const std::string a_outputRealm, const int a_level, const bool a_interpToCentroids, const bool a_interpGhost) const noexcept
 Write data to output. Convenience function. More...
 

Protected Attributes

bool m_transparentEB
 Turn on/off transparent boundaries.
 
bool m_instantaneous
 Instantaneous transport or not.
 
bool m_blendConservation
 Flag for blending the deposition clouds with the nonconservative divergence.
 
bool m_depositNumber
 If true, the NUMBER of of Photons will be deposited in each cell.
 
bool m_plotNumbers
 Switch for plotting numbers or densities.
 
bool m_plotPhotons
 Check if m_photons should be plotted.
 
bool m_plotBulkPhotons
 Check if m_bulkPhotons should be plotted.
 
bool m_plotSourcePhotons
 Check if source_bulkPhotons should be plotted.
 
bool m_plotEBPhotons
 Check if m_ebPhotons should be plotted.
 
bool m_plotDomainPhotons
 Check if m_ebPhotons should be plotted.
 
bool m_dirtySampling
 Dirty sampling or not.
 
size_t m_maxPhotonsGeneratedPerCell
 Number of computational photons generated per cell.
 
int m_numSamplingPackets
 
int m_seed
 RNG seed.
 
int m_haloBuffer
 Halo size for particles.
 
Real m_bisectStep
 Bisection step for when we use the bisection intersection algorithm.
 
PhotonGeneration m_photoGenerationMethod
 Photon generation type.
 
SourceType m_sourceType
 Source type.
 
IntersectionEB m_intersectionEB
 
DepositionType m_deposition
 Deposition type.
 
CoarseFineDeposition m_coarseFineDeposition
 Coarse-fine deposition strategy.
 
EBAMRCellData m_scratch
 Coarse data for interpolation of deposition clouds.
 
EBAMRIVData m_depositionNC
 Scratch storage for holding the non-conservative deposition.
 
EBAMRIVData m_massDiff
 Scratch storage for holding the mass difference when using hybrid deposition.
 
ParticleContainer< Photonm_photons
 All particles.
 
ParticleContainer< Photonm_bulkPhotons
 Photons absorbed in the volume.
 
ParticleContainer< Photonm_ebPhotons
 This is a particle container for Photons that crossed EBs. It is filled during the advance step.
 
ParticleContainer< Photonm_domainPhotons
 This is a particle container for Photons that crossed boundaries. It is filled during the advance step.
 
ParticleContainer< Photonm_sourcePhotons
 This is a particle container that can be used to add Photons directly.
 
- Protected Attributes inherited from RtSolver
Location::Cell m_dataLocation
 Data location.
 
std::string m_realm
 Realm where this solver lives.
 
RefCountedPtr< EBIndexSpace > m_ebis
 EBIndexSpace for this solver.
 
RefCountedPtr< RtSpeciesm_rtSpecies
 Radiative transfer species (contains meta-information like initial conditions)
 
RefCountedPtr< ComputationalGeometrym_computationalGeometry
 Computational geometry.
 
RefCountedPtr< AmrMeshm_amr
 AMR; needed for grid stuff.
 
phase::which_phase m_phase
 Phase.
 
std::string m_name
 Name for this solver.
 
std::string m_className
 Class name – needed because inherited classes will be named different.
 
EBAMRCellData m_cachePhi
 Cached state used for regridding.
 
EBAMRCellData m_phi
 Internal state. More...
 
EBAMRCellData m_source
 Source term. More...
 
EBAMRFluxData m_kappa
 Absorption coefficient.
 
EBAMRIVData m_kappaEB
 Absorption coefficient on EB faces.
 
Real m_time
 Time.
 
Real m_dt
 Time increment.
 
bool m_stationary
 Stationary solver or not.
 
bool m_plotPhi
 Output state.
 
bool m_plotSource
 Output source term.
 
int m_verbosity
 
int m_timeStep
 Time step.
 

Additional Inherited Members

- Static Protected Attributes inherited from RtSolver
static constexpr int m_comp = 0
 Default component that we solve for.
 
static constexpr int m_nComp = 1
 Default number of components.
 

Detailed Description

Radiative tranfer equation solver using Monte-Carlo simulation.

This class is, by default, a non-stationary radiative transfer solver which simulates computational photons that can also be deposited on the mesh. Various options are available for specifying how the photons are generated, and how they are propagated (e.g. with speed of light or instantaneously). Note that this solver runs its own intersection tests with the EB due to dependencies on how the photons are moved (we store free-flight photons, absorbed photons, intersected photons, etc).

Note
The solver is missing appropriate boundary conditions on domain faces.

Member Enumeration Documentation

◆ IntersectionEB

enum McPhoto::IntersectionEB
strongprotected

An enum for switching between various types of EB intersection algorithms when intersecting photons with the EB.

Raycast means ray-casting algorithm. Bisection means that the traveled path is divided into intervals and we apply a bisection algorithm for computing the intersection point. LSF is just ray-casting, but it uses the implicit function on the mesh rather than calling it directly.

Member Function Documentation

◆ advance()

bool McPhoto::advance ( const Real  a_dt,
EBAMRCellData a_phi,
const EBAMRCellData a_source,
const bool  a_zeroPhi = false 
)
overridevirtual

Advance RTE and deposit photon particles on the mesh.

Parameters
[in]a_dtTime step
[out]a_phiPhoton density on the mesh
[in]a_sourceSource term, i.e. number of photons produced per unit time. The photons will be generated based on what the user specifies in the input parameters for this class. I.e., whether or not a_source will generate stochastic or deterministic photons, and scaled based on what a_source actually contains (number of photons, number/volume, number/second, number/(volume*second)).
[in]a_zeroPhiDead parameter not used in this implementation.
Note
This is the "fluid interface" to the radiative transfer solver in which the solver uses discrete photons for sampling the isotropic density. Particle codes will probably use advancePhotonsInstantaneous or advancePhotonsTransient .

This routine switches between stationary and transient solvers depending on whether or not the solver is stationary.

Implements RtSolver.

◆ advancePhotonsInstantaneous()

void McPhoto::advancePhotonsInstantaneous ( ParticleContainer< Photon > &  a_bulkPhotons,
ParticleContainer< Photon > &  a_ebPhotons,
ParticleContainer< Photon > &  a_domainPhotons,
ParticleContainer< Photon > &  a_photons 
)
virtual

Move photons and absorb them on various objects.

Parameters
[out]a_bulkPhotonsPhotons absorbed on the mesh
[out]a_ebPhotonsPhotons absorbed on the EB
[out]a_domainPhotonsPhotons absorbed on the domain edges (faces)
[in,out]a_photonsOriginal photons
Note
This routine moves photons with an instantaneous kernel, i.e. all photons are always absorbed and none are left behind as free-flight photons (i.e. a_photons will be empty on output)

◆ advancePhotonsTransient()

void McPhoto::advancePhotonsTransient ( ParticleContainer< Photon > &  a_bulkPhotons,
ParticleContainer< Photon > &  a_ebPhotons,
ParticleContainer< Photon > &  a_domainPhotons,
ParticleContainer< Photon > &  a_photons,
const Real  a_dt 
)
virtual

Move photons and absorb them on various objects.

Parameters
[out]a_bulkPhotonsPhotons absorbed on the mesh
[out]a_ebPhotonsPhotons absorbed on the EB
[out]a_domainPhotonsPhotons absorbed on the domain edges (faces)
[in,out]a_photonsOriginal photons
[in]a_dtTime step
Note
This routine moves photons with a transient kernel, i.e. all photons propagate at most c*dt and checks are made to determine if they are absorbed on the mesh or various objects (EB, domain). a_photons may be non-empty on output, in which case there are still free-flight photons in the solver.

◆ clear() [1/3]

void McPhoto::clear ( AMRParticles< Photon > &  a_photons)
virtual

Clear data holder.

Parameters
[in]a_photonsWhich container to clear

◆ clear() [2/3]

void McPhoto::clear ( const WhichContainer a_which)
virtual

Clear data holder.

Parameters
[in]a_whichWhich container to clear

◆ clear() [3/3]

void McPhoto::clear ( ParticleContainer< Photon > &  a_photons)
virtual

Clear data holder.

Parameters
[in]a_photonsWhich container to clear

◆ computeBoundaryFlux()

void McPhoto::computeBoundaryFlux ( EBAMRIVData a_ebFlux,
const EBAMRCellData a_phi 
)
overridevirtual

Compute the boundary flux.

Parameters
[out]a_ebFluxFlux on the EB
[in]a_phiRTE solution
Note
This sets the boundary flux to zero. If you want to use a boundary flux, use it through m_ebPhotons

Implements RtSolver.

◆ computeDensity()

void McPhoto::computeDensity ( EBAMRCellData a_isotropic,
const EBAMRCellData a_phi 
)
overridevirtual

Compute isotropic radiative density from mesh solution.

Parameters
[out]a_isotropicIsotropic part of RTE solution
[in]a_phiRTE solution
Note
Issues an error - this routine doesn't make a lot of sense, really.

Implements RtSolver.

◆ computeDomainFlux()

void McPhoto::computeDomainFlux ( EBAMRIFData a_domainFlux,
const EBAMRCellData a_phi 
)
overridevirtual

Compute the domain flux.

Parameters
[out]a_domainFluxDomain flux
[in]a_phiRTE solution
Note
This sets the boundary flux to zero. If you want to use a domain flux, use it through m_domainPhotons

Implements RtSolver.

◆ computeFlux()

void McPhoto::computeFlux ( EBAMRCellData a_flux,
const EBAMRCellData a_phi 
)
overridevirtual

Compute the flux.

Parameters
[out]a_fluxFlux
[in]a_phiRTE solution
Note
Issues an error – currently don't know how to compute the mesh flux from computational photons....

Implements RtSolver.

◆ computeLoads()

void McPhoto::computeLoads ( Vector< long long > &  a_loads,
const DisjointBoxLayout &  a_dbl,
const int  a_level 
) const
overridevirtualnoexcept

Get computational loads for a specific grid level. This computes the number of photons in the bulk data holder.

Parameters
[out]a_loadsGrid loads for this level.
[in]a_dblGrids on input level
[in]a_levelInput level
Returns
Loads for each box on a grid level.
Note
The return vector should have the same order as the boxes in a_dbl. E.g. ret[0] must be the load for a_dbl.boxArray()[0];

The default implementation returns the number of cells in the grid patch as a proxy for the load.

Reimplemented from RtSolver.

◆ computeNumPhysicalPhotons()

void McPhoto::computeNumPhysicalPhotons ( EBAMRCellData a_numPhysPhotonsTotal,
EBAMRCellData a_numPhysPhotonsPerPacket,
const EBAMRCellData a_source,
const Real  a_dt 
) const
virtualnoexcept

Compute the number of physical photons in each grid cell.

Parameters
[out]a_numPhysPhotonsTotalTotal number of physical photons to be generated in each grid cell.
[out]a_numPhysPhotonsPacketNumber of physical photons to be generated in each packet. Must have a components = m_numSamplingPackets
[in]a_sourceSource term – the interpretation of this depends on the user-specified interpretation of the source term.
[in]a_dtTime step

◆ countPhotons()

int McPhoto::countPhotons ( const AMRParticles< Photon > &  a_photons) const
virtual

Count number of photons in particle list.

Parameters
[in]a_photonsParticle list

◆ depositHybrid()

void McPhoto::depositHybrid ( EBAMRCellData a_depositionH,
EBAMRIVData a_massDifference,
const EBAMRIVData a_depositionNC 
) const
protectednoexcept

Make the hybrid deposition. Also compute the mass difference.

Parameters
[in,out]a_depositionHOn input this should be the conservative deposition
[out]a_massDifferenceMass difference between conservative and non-conservative
[in]a_depositioNCNon-conservative deposition

◆ depositKappaConservative()

template<class P , const Real &(P::*)() const particleScalarField>
void McPhoto::depositKappaConservative ( EBAMRCellData a_phi,
ParticleContainer< P > &  a_particles,
const DepositionType  a_deposition,
const CoarseFineDeposition  a_coarseFineDeposition 
) const
protectednoexcept

This computes the "conservative" deposition, multiplied by kappa.

Parameters
[out]a_phiMesh density.
[in]a_particlesParticles to be deposited
[in]a_depositionDeposition method
[in]a_coarseFineDepositionCoarse-fine deposition strategy.

◆ depositNonConservative()

void McPhoto::depositNonConservative ( EBAMRIVData a_depositionNC,
const EBAMRCellData a_depositionKappaC 
) const
protectednoexcept

Make the "non-conservative" kappa deposition.

Parameters
[out]a_depositionNCNon-conservative deposition in cut-cells
[in]a_depositionKappaCConservative deposition

◆ depositPhotons() [1/2]

void McPhoto::depositPhotons ( )
virtual

Deposit photons on the mesh.

This deposits m_photons onto m_phi

◆ depositPhotons() [2/2]

template<class P , const Real &(P::*)() const particleScalarField>
void McPhoto::depositPhotons ( EBAMRCellData a_phi,
ParticleContainer< P > &  a_particles,
const DepositionType a_deposition 
) const
noexcept

Deposit photons.

Parameters
[out]a_phiDeposited mesh density
[in,out]a_photonsComputational photons to be deposited
[in]a_depositionDeposition method

◆ depositPhotonsNGP()

void McPhoto::depositPhotonsNGP ( LevelData< EBCellFAB > &  a_output,
const ParticleContainer< Photon > &  a_particles,
const int  a_level 
) const
protectedvirtualnoexcept

Do an NGP deposit on a specific grid level. Used for IO.

Parameters
[out]a_outputContains NGP deposition of photons. Ignores cut-cells.
[in]a_photonsphotons
[in]a_levelGrid level

◆ dirtySamplePhotons()

void McPhoto::dirtySamplePhotons ( ParticleContainer< PointParticle > &  a_photons,
EBAMRCellData a_phi,
const EBAMRCellData a_numPhysicalPhotons,
const size_t  a_maxPhotonsPerCell 
) const
virtualnoexcept

Dirty-sampling method for photons.

This is a "dirty" routine in the sense that it does not respect the EB and avoids filling m_bulkParticles. If you also use an NGP scheme this avoids communication between sampling buckets. But this is not something that we recommend using, unless you have a geometry-less domain. If those are good for you, feel free to use this routine.

Parameters
[in,out]a_photonsComputational photons. Will be populated and depleted inside the routine.
[in,out]a_phiMesh-based density. Will be incremented by the generated photons.
[in]a_numPhysicalPhotonsNumber of physical photons to add to a_photons
[in]a_maxPhotonsPerCellMaximum number of photons generated per cell.

◆ domainBcMap()

int McPhoto::domainBcMap ( const int  a_dir,
const Side::LoHiSide  a_side 
)
protected

Mapping function for domain boundary conditions.

Parameters
[in]a_dirCoordinate direction
[in]a_sideCoordinate side

◆ drawPhotons()

size_t McPhoto::drawPhotons ( const Real  a_source,
const Real  a_volume,
const Real  a_dt 
) const
protectednoexcept

Draw photons in a cell and volume.

Parameters
[in]a_sourceSource term, i.e. number of photons generated/x
[in]a_volumeGrid cell volume
[in]a_dtTime step

The return result of this function depends on the input parameters to this class.

◆ generateComputationalPhotons()

void McPhoto::generateComputationalPhotons ( ParticleContainer< Photon > &  a_photons,
const EBAMRCellData a_numPhysicalPhotons,
const size_t  a_maxPhotonsPerCell 
) const
virtualnoexcept

Generate computational photons.

Parameters
[in,out]a_photonsComputational photons
[in]a_numPhysicalPhotonsNumber of physical photons to add to a_photons
[in]a_maxPhotonsPerCellMaximum number of photons generated per cell.

◆ getBulkPhotons()

ParticleContainer< Photon > & McPhoto::getBulkPhotons ( )
virtual

Get bulk photons, i.e. photons absorbed on the mesh.

Returns
m_bulkPhotons

◆ getDomainPhotons()

ParticleContainer< Photon > & McPhoto::getDomainPhotons ( )
virtual

Get domain photons, i.e. photons absorbed on domain edges/faces.

Returns
m_domainPhotonsl

◆ getEbPhotons()

ParticleContainer< Photon > & McPhoto::getEbPhotons ( )
virtual

Get eb Photons, i.e. photons absorbed on the EB.

Returns
m_ebPhotons

◆ getNumberOfPlotVariables()

int McPhoto::getNumberOfPlotVariables ( ) const
overridevirtual

Get number of output variables.

Returns
Returns number of plot variables that will be written to file in writePlotData

Reimplemented from RtSolver.

◆ getPhotons()

ParticleContainer< Photon > & McPhoto::getPhotons ( )
virtual

Get m_photons.

Returns
m_photons

◆ getPlotVariableNames()

Vector< std::string > McPhoto::getPlotVariableNames ( ) const
overridevirtual

Write checkpoint data into handle.

Parameters
[in,out]a_handleHDF5 file
[in]a_levelGrid level

Read checkpoint data from handle

Parameters
[in,out]a_handleHDF5 file
[in]a_levelGrid level

Get output plot names

Returns
Returns a list of plot variable names.

Reimplemented from RtSolver.

◆ getSourcePhotons()

ParticleContainer< Photon > & McPhoto::getSourcePhotons ( )
virtual

Get source photons.

Returns
m_sourcePhotons

◆ parseDirtySampling()

void McPhoto::parseDirtySampling ( )
protected

Parse dirty sampling for photons.

This is a hidden option – I really don't want it to be part of the regular user interface because the code is not general (only works for instantaneous, non-EB cases where photons aren't even absorbed on the domain boundaries)

◆ preRegrid()

void McPhoto::preRegrid ( const int  a_base,
const int  a_oldFinestLevel 
)
overridevirtual

preRegrid operations

Parameters
[in]a_baseCoarsest level which changes during regrid
[in]a_oldFinestLevelFinest grid level before regrid

Implements RtSolver.

◆ randomExponential()

Real McPhoto::randomExponential ( const Real  a_mean) const
protectednoexcept

Random exponential trial.

Parameters
[in]a_meanMean value in exponential distribution.

◆ regrid()

void McPhoto::regrid ( const int  a_lmin,
const int  a_oldFinestLevel,
const int  a_newFinestLevel 
)
overridevirtual

Regrid function for this class.

Parameters
[in]a_lminCoarsest level that changed during regrid
[i]a_oldFinestLevel Finest grid level before the regrid
[i]a_newFinestLevel Finest grid level after the regrid

Implements RtSolver.

◆ remap()

void McPhoto::remap ( ParticleContainer< Photon > &  a_photons)
virtual

Remap computational particles.

Parameters
[in]a_photonsComputational particles to be remapped.

◆ sortPhotonsByCell()

void McPhoto::sortPhotonsByCell ( const WhichContainer a_which)
virtual

Sort container by cell.

WhichContainer::Photon = m_photons, WhichContainer::EB = m_ebPhotons and so on.

Parameters
[in]a_whichWhich container to sort.

◆ sortPhotonsByPatch()

void McPhoto::sortPhotonsByPatch ( const WhichContainer a_which)
virtual

Sort container by patch.

WhichContainer::Photon = m_photons, WhichContainer::EB = m_ebPhotons and so on.

Parameters
[in]a_whichWhich container to sort.

◆ writePlotData()

void McPhoto::writePlotData ( LevelData< EBCellFAB > &  a_output,
int &  a_icomp,
const std::string  a_outputRealm,
const int  a_level 
) const
overridevirtualnoexcept

Write plot data.

Parameters
[in,out]a_ouputOutput data holder
[in,out]a_icompStarting component where we begin to write data to a_output
[in]a_outputRealmRealm where a_output belongs
[in]a_levelGrid level

Reimplemented from RtSolver.

◆ writePlotFile()

void McPhoto::writePlotFile ( )
overridevirtual

Write plot file.

Note
Currently calls an error because it is not implemented.

Implements RtSolver.


The documentation for this class was generated from the following files: