ABLATE Source Documentation  0.12.34
ablate::radiation::OrthogonalRadiation Class Reference
+ Inheritance diagram for ablate::radiation::OrthogonalRadiation:

Public Member Functions

 OrthogonalRadiation (const std::string &solverId, const std::shared_ptr< domain::Region > &region, std::shared_ptr< eos::radiationProperties::RadiationModel > radiationModelIn, std::shared_ptr< ablate::monitors::logs::Log >={})
 
 ~OrthogonalRadiation ()
 The ray number should never be used because there is only one ray emanating from every boundary face.
 
void Setup (const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain) override
 
void GetSurfaceIntensity (PetscReal *intensityReturn, PetscInt faceId, PetscReal temperature, PetscReal emissivity=1.0, PetscReal absorptivity=1.0) override
 
- Public Member Functions inherited from ablate::radiation::SurfaceRadiation
 SurfaceRadiation (const std::string &solverId, const std::shared_ptr< domain::Region > &region, const PetscInt raynumber, std::shared_ptr< eos::radiationProperties::RadiationModel > radiationModelIn, std::shared_ptr< ablate::monitors::logs::Log >={})
 
void Initialize (const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain) override
 
PetscReal SurfaceComponent (const PetscReal normal[], PetscInt iCell, PetscInt nphi, PetscInt ntheta) override
 
- Public Member Functions inherited from ablate::radiation::Radiation
 Radiation (const std::string &solverId, const std::shared_ptr< domain::Region > &region, const PetscInt raynumber, std::shared_ptr< eos::radiationProperties::RadiationModel > radiationModelIn, std::shared_ptr< ablate::monitors::logs::Log >={})
 
void GetIntensity (PetscReal *intensity, PetscInt index, const domain::Range &cellRange, PetscReal temperature, PetscReal kappa)
 
std::string GetId ()
 
eos::ThermodynamicTemperatureFunction GetAbsorptionFunction ()
 
void EvaluateGains (Vec solVec, ablate::domain::Field temperatureField, Vec auxVec)
 
virtual void ParticleStep (ablate::domain::SubDomain &subDomain, DM faceDM, const PetscScalar *faceGeomArray, DM radReturn, PetscInt nlocalpoints, PetscInt nglobalpoints)
 Routine to move the particle one step. More...
 
virtual void IdentifyNewRaysOnRank (domain::SubDomain &subDomain, DM radReturn, PetscInt nlocalpoints)
 If this local rank has never seen this search particle before, then it needs to add a new ray segment to local memory and record its index. More...
 
std::shared_ptr< eos::radiationProperties::RadiationModelGetRadiationModel ()
 provide access to the model used to provided the absorptivity function
 
std::shared_ptr< domain::RegionGetRegion () const
 

Static Public Member Functions

static std::string GetClassType ()
 
- Static Public Member Functions inherited from ablate::radiation::SurfaceRadiation
static std::string GetClassType ()
 
- Static Public Member Functions inherited from ablate::radiation::Radiation
static PetscReal GetBlackBodyTotalIntensity (PetscReal temperature, const PetscReal refractiveIndex)
 
static PetscReal GetBlackBodyWavelengthIntensity (const PetscReal temperature, const PetscReal wavelength, const PetscReal refractiveIndex)
 
static std::string GetClassType ()
 

Additional Inherited Members

- Protected Member Functions inherited from ablate::radiation::Radiation
PetscReal FaceIntersect (PetscInt ip, Virtualcoord *virtualcoord, PetscFVFaceGeom *face) const
 Class Methods. More...
 
void UpdateCoordinates (PetscInt ipart, Virtualcoord *virtualcoord, PetscReal *coord, PetscReal adv) const
 
void DeleteOutOfBounds (ablate::domain::SubDomain &subDomain)
 
virtual void SetBoundary (CellSegment &raySegment, PetscInt index, Identifier identifier)
 
- Protected Member Functions inherited from ablate::utilities::Loggable< Radiation >
const PetscClassId & GetPetscClassId () const
 
PetscLogEvent RegisterEvent (const char *eventName)
 
void StartEvent (const char *eventName) const
 
void EndEvent () const
 
- Protected Attributes inherited from ablate::radiation::SurfaceRadiation
ablate::domain::ReverseRange indexLookup
 used to look up from the face id to range index
 
- Protected Attributes inherited from ablate::radiation::Radiation
DM radSearch = nullptr
 DM which the search particles occupy. This representations the physical particle in space.
 
Vec faceGeomVec = nullptr
 Vector used to describe the entire face geom of the dm. This is constant and does not depend upon region.
 
Vec cellGeomVec = nullptr
 
MPI_Datatype carrierMpiType
 create a data type to simplify moving the carrier
 
PetscInt dim = 0
 Class inputs and Variables. More...
 
PetscInt nTheta
 The number of angles to solve with, given by user input.
 
PetscInt nPhi
 The number of angles to solve with, given by user input (x2)
 
PetscReal minCellRadius {}
 
std::vector< std::vector< CellSegment > > raySegments
 store the local rays identified on this rank. This includes rays that do and do not originate on this rank
 
std::vector< CarrierraySegmentsCalculations
 the calculation over each of the remoteRays. indexed over remote ray
 
PetscInt numberOriginRays
 store the number of originating rays
 
PetscInt numberOriginCells
 store the number of originating rays cells
 
PetscInt raysPerCell
 the number of rays per cell
 
std::vector< unsigned short int > raySegmentsPerOriginRay
 store the number of ray segments for each originating on this rank. This may be zero
 
std::vector< CarrierraySegmentSummary
 a vector of raySegment information for every local/remote ray segment ordered as ray, segment
 
std::vector< PetscReal > gainsFactor
 the factor for each origin ray
 
std::vector< PetscScalar > evaluatedGains
 size up the evaluated gains, this index is based upon order of the requested cells
 
PetscSF remoteAccess = nullptr
 Store the petscSF that is used for pulling remote ray calculation.
 
std::string solverId
 the name of this solver
 
const std::shared_ptr< domain::Regionregion
 the region for which this solver applies
 
const std::shared_ptr< eos::radiationProperties::RadiationModelradiationModel
 model used to provided the absorptivity function
 
eos::ThermodynamicTemperatureFunction absorptivityFunction
 hold a pointer to the absorptivity function
 
eos::ThermodynamicTemperatureFunction emissivityFunction
 
const std::shared_ptr< ablate::monitors::logs::Loglog = nullptr
 
- Static Protected Attributes inherited from ablate::radiation::Radiation
static constexpr char IdentifierField [] = "identifier"
 
static constexpr char VirtualCoordField [] = "virtual coord"
 

Member Function Documentation

◆ GetClassType()

static std::string ablate::radiation::OrthogonalRadiation::GetClassType ( )
inlinestatic

Represents the name of the class for logging and other utilities

Returns

◆ GetSurfaceIntensity()

void ablate::radiation::OrthogonalRadiation::GetSurfaceIntensity ( PetscReal *  intensity,
PetscInt  faceId,
PetscReal  temperature,
PetscReal  emissivity = 1.0,
PetscReal  absorptivity = 1.0 
)
inlineoverridevirtual

Compute total intensity (pre computed gains + current loss) with

Parameters
faceIdthe current face id
temperaturethe temperature of the face
emissivitythe emissivity of the surface
Returns

Reimplemented from ablate::radiation::SurfaceRadiation.

◆ Setup()

void ablate::radiation::OrthogonalRadiation::Setup ( const ablate::domain::Range cellRange,
ablate::domain::SubDomain subDomain 
)
overridevirtual

SubDomain Register and Setup

allows initialization after the subdomain and dm is established

< Number of dimensions already defined in the setup

< Reduce the number of rays if one dimensional symmetry can be taken advantage of

Begins radiation properties model Runs the ray initialization, finding cell indices Initialize the log if provided

Initialization to call, draws each ray vector and gets all of the cells associated with it (sorted by distance and starting at the boundary working in) This is done by creating particles at the center of each cell and iterating through them Get setup things for the position vector of the current cell index Declare the variables that will contain the geometry of the cells Obtain the geometric information about the cells in the DM

do a simple sanity check for labels

< Get the origin rank of the current process. The particle belongs to this rank. The rank only needs to be read once.

Setup the particles and their associated fields including: origin domain/ ray identifier / # domains crossed, and coordinates. Instantiate ray particles for each local cell only.

< Number of points to insert into the particle field. One particle for each ray.

Create the DMSwarm

Configure radsearch to be of type PIC/Basic

Register fields within the DMSwarm

< A field to store the ray identifier [origin][iCell][ntheta][nphi][ndomain]

< A field representing the three dimensional coordinates of the particle. Three "virtual" dims are required.

< Initialize the fields that have been defined

Set initial local sizes of the DMSwarm with a buffer length of zero

< Set the number of initial particles to the number of rays in the subdomain. Set the buffer size to zero.

Declare some information associated with the field declarations

< Pointer to the coordinate field information

< Pointer to the cell index information

< Pointer to the primary (virtual) coordinate field information

< Pointer to the ray identifier information

Get the fields associated with the particle swarm so that they can be modified

< Initialize a counter to represent the particle index. This will be iterated every time that the inner loop is passed through.

< This will iterate only though local cells

< Isolates the valid cells

Since we don't know which side of the face is inside the region, we can just spawn two particles on opposite sides of the face and give them both the same gains factor weighting. One of them will be on the correct side and one if them will not. Whichever one is not will be deleted. The iteration represents a loop where i is the coefficient of the normal vector. It is equal to -1 on the first iteration and equal to 1 on the second iteration. This will put a particle on both sides of the normal vector.

Update the direction vector of the search particle

< x component direction from the cell face normal

< y component direction from the cell face normal

< z component direction from the cell face normal

Get the particle coordinate field and write the cellGeom->centroid[xyz] into it

< Offset from the centroid slightly so they sit in a cell if they are on its face.

Update the physical coordinate field so that the real particle location can be updated.

Update the physical coordinate field so that the real particle location can be updated.

adv value of 0.0 places the particle exactly where the virtual coordinates are.

Label the particle with the ray identifier. With what is known at this point

< Input the ray identifier. This location scheme represents stepping four entries for every particle index increase

This serve to identify the ray id on the origin

Intensity of the irradiation at the surface is pi * intensity integrated.

Set the index of the field value so that it can be written to for every particle

< Must be iterated at the end since the value is initialized at zero.}

Restore the fields associated with the particles

< Sets the search particles in the cell indexes to which they have been assigned

< Number of dimensions already defined in the setup

< Reduce the number of rays if one dimensional symmetry can be taken advantage of

Begins radiation properties model Runs the ray initialization, finding cell indices Initialize the log if provided

Initialization to call, draws each ray vector and gets all of the cells associated with it (sorted by distance and starting at the boundary working in) This is done by creating particles at the center of each cell and iterating through them Get setup things for the position vector of the current cell index Declare the variables that will contain the geometry of the cells Obtain the geometric information about the cells in the DM

do a simple sanity check for labels

< Get the origin rank of the current process. The particle belongs to this rank. The rank only needs to be read once.

Setup the particles and their associated fields including: origin domain/ ray identifier / # domains crossed, and coordinates. Instantiate ray particles for each local cell only.

< Number of points to insert into the particle field. One particle for each ray.

Create the DMSwarm

Configure radsearch to be of type PIC/Basic

Register fields within the DMSwarm

< A field to store the ray identifier [origin][iCell][ntheta][nphi][ndomain]

< A field representing the three dimensional coordinates of the particle. Three "virtual" dims are required.

< Initialize the fields that have been defined

Set initial local sizes of the DMSwarm with a buffer length of zero

< Set the number of initial particles to the number of rays in the subdomain. Set the buffer size to zero.

Declare some information associated with the field declarations

< Pointer to the coordinate field information

< Pointer to the cell index information

< Pointer to the primary (virtual) coordinate field information

< Pointer to the ray identifier information

Get the fields associated with the particle swarm so that they can be modified

< Initialize a counter to represent the particle index. This will be iterated every time that the inner loop is passed through.

< This will iterate only though local cells

< Isolates the valid cells

for every angle theta for every angle phi

Get the initial direction of the search particle from the angle number that it was initialized with

< Theta angle of the ray

< Phi angle of the ray

Update the direction vector of the search particle

< x component conversion from spherical coordinates, adding the position of the current cell

< y component conversion from spherical coordinates, adding the position of the current cell

< z component conversion from spherical coordinates, adding the position of the current cell

Get the particle coordinate field and write the cellGeom->centroid[xyz] into it

< Offset from the centroid slightly so they sit in a cell if they are on its face.

Update the physical coordinate field so that the real particle location can be updated.

adv value of 0.0 places the particle exactly where the virtual coordinates are.

Label the particle with the ray identifier. With what is known at this point

< Input the ray identifier. This location scheme represents stepping four entries for every particle index increase

This serve to identify the ray id on the origin

Set the index of the field value so that it can be written to for every particle

< Must be iterated at the end since the value is initialized at zero.

Restore the fields associated with the particles

< Sets the search particles in the cell indexes to which they have been assigned

Reimplemented from ablate::radiation::Radiation.


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