ABLATE Source Documentation  0.12.33
ablate::radiation::Radiation Class Reference
+ Inheritance diagram for ablate::radiation::Radiation:

Classes

struct  Carrier
 
struct  CellSegment
 
struct  Identifier
 
struct  Virtualcoord
 

Public Member Functions

 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 >={})
 
virtual void Setup (const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain)
 
virtual void Initialize (const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain)
 
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...
 
virtual PetscReal SurfaceComponent (const PetscReal normal[], PetscInt iCell, PetscInt nphi, PetscInt ntheta)
 
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 PetscReal GetBlackBodyTotalIntensity (PetscReal temperature, const PetscReal refractiveIndex)
 
static PetscReal GetBlackBodyWavelengthIntensity (const PetscReal temperature, const PetscReal wavelength, const PetscReal refractiveIndex)
 
static std::string GetClassType ()
 

Protected Member Functions

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

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

static constexpr char IdentifierField [] = "identifier"
 
static constexpr char VirtualCoordField [] = "virtual coord"
 

Constructor & Destructor Documentation

◆ Radiation()

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 log = {} 
)

< Cell solver provides cell based functionality, right hand side function compatibility with finite element/ volume, loggable allows for the timing and tracking of events

Parameters
solverIdthe id for this solver
regionthe boundary cell region
rayNumber
optionsother options

Member Function Documentation

◆ DeleteOutOfBounds()

void ablate::radiation::Radiation::DeleteOutOfBounds ( ablate::domain::SubDomain subDomain)
protected

Delete the particles that are outside of the domain This is used in the initialization for ray tracing solvers which rely on domain bounds to inform the placement of particles around the boundary regions.

Parameters
subDomain

< Pointer to the coordinate field 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

< Recalculate the number of particles that are in the domain

< If the particles that were just created are sitting in the boundary cell of the face that they belong to, delete them

< If the particle location index and boundary cell index are the same, then they should be deleted

< Delete the particle!

< Check the point replacing the one that was deleted

Restore the fields associated with the particles

◆ EvaluateGains()

void ablate::radiation::Radiation::EvaluateGains ( Vec  solVec,
ablate::domain::Field  temperatureField,
Vec  auxVec 
)

Evaluates the ray intensity from the domain to update the effects of irradiation. Does not impact the solution unless the solve function is called again.

Get the array of the solution vector.

Get the array of the aux vector.

Zero this ray segment for all wavelengths

Iterate through every wavelength entry in this ray segment

This is allowed to be cast to auto and indexed raySegmentIndex because there is only one physical ray segment that we are reading from.

< The solution value at any given location

< The temperature at any given location

Input absorptivity (kappa) values from model here.

< Absorptivity coefficient, property of each cell. This is an array that we will iterate through for every evaluation

Get the absorption and emission information from the provided properties models.

Get the pointer to the returned array of absorption values. Iterate through every wavelength for the evaluation.

In the future we may want to set this intensity with a boundary condition class.

March over each INDEXING ANNOTATIONS: evaluatedGains: There is a gain evaluation for every cell * wavelength Therefore, the indexing of the gain evaluation is [absorptivityFunction.propertySize * cell + i] because the cells are looped through outside of the wavelengths raySegmentsPerOriginRay: This only stores the number of segments in each ray. There is no reason to index this with wavelength. Therefore, the indexing is [rayOffset], where the ray refers to the wavelength independent ray count. raySegmentSummary: This will store a value for every ray segment and wavelength. Each ray will integrate its ray segments together for every wavelength.

Zero the evaluated gains for this ray specifically. Do this for all wavelengths.

for each segment in this ray Integrate the wavelength dependent intensity calculation for each segment for each wavelength We need to store all of the wavelength results on the final evaluation of the cell Therefore, we should first iterate through the wavelengths first and sum the effects of every wavelength on every cell.

Cleanup

◆ FaceIntersect()

PetscReal ablate::radiation::Radiation::FaceIntersect ( PetscInt  ip,
Virtualcoord virtualcoord,
PetscFVFaceGeom *  face 
) const
protected

Class Methods.

Returns the forward path length of a travelling particle with any face. The function will return zero if the intersection is not in the direction of travel. x

Parameters
virtualcoordthe struct containing particle position information
facethe struct containing information about a cell face Returns the distance away from a virtual coordinate at which its path intersects a line.

<(planeNormal.dot(planePoint) - planeNormal.dot(linePoint)) / planeNormal.dot(lineDirection.normalize())

◆ GetBlackBodyTotalIntensity()

static PetscReal ablate::radiation::Radiation::GetBlackBodyTotalIntensity ( PetscReal  temperature,
const PetscReal  refractiveIndex 
)
inlinestatic

Returns the black body intensity for a given temperature and emissivity

Parameters
temperature
refractiveIndex
Returns

◆ GetBlackBodyWavelengthIntensity()

static PetscReal ablate::radiation::Radiation::GetBlackBodyWavelengthIntensity ( const PetscReal  temperature,
const PetscReal  wavelength,
const PetscReal  refractiveIndex 
)
inlinestatic

Returns the black body intensity at a specific wavelength. This must be integrated over a bandwidth.

Parameters
temperature
wavelength
refractiveIndex
Returns

◆ GetClassType()

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

Represents the name of the class for logging and other utilities

Returns

◆ GetIntensity()

void ablate::radiation::Radiation::GetIntensity ( PetscReal *  intensity,
PetscInt  index,
const domain::Range cellRange,
PetscReal  temperature,
PetscReal  kappa 
)
inline

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

Parameters
indexthe current index in the cell range (note, this is not the cell/face id)
cellRangethe cell/face range used in setup/initialize
temperaturethe temperature of the cell or face
kappathe absorptivity of the cell
Returns

This should take the emission as an input (if non-grey)

◆ GetRegion()

std::shared_ptr<domain::Region> ablate::radiation::Radiation::GetRegion ( ) const
inline

public function to return the region for this radiation solver

Returns

◆ IdentifyNewRaysOnRank()

void ablate::radiation::Radiation::IdentifyNewRaysOnRank ( domain::SubDomain subDomain,
DM  radReturn,
PetscInt  nlocalpoints 
)
virtual

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.

Check that the particle is in a valid region

Declare some information associated with the field declarations

< Pointer to the ray identifier information

Get all of the ray information from the particle Get the ntheta and nphi from the particle that is currently being looked at. This will be used to identify its ray and calculate its direction.

Check that the particle is in a valid region

< Another solve particle is added here because the search particle has entered a new domain

< Pointer to the ray identifier information

while we are here, set the return rank. This won't change anything until migrate is called

< Get the fields from the radsolve swarm so the new point can be written to them

< Get the fields from the radsolve swarm so the new point can be written to them

Reimplemented in ablate::radiation::RaySharingRadiation.

◆ Initialize()

void ablate::radiation::Radiation::Initialize ( const ablate::domain::Range cellRange,
ablate::domain::SubDomain subDomain 
)
virtual
Parameters
cellRangeThe range of cells for which rays are initialized

< Initialize the fields that have been defined

This will be added to as rays are created on each rank

< Get the geometry vectors

Exact some information associated with the field declarations from the swarm

< Pointer to the coordinate field information

< Pointer to the primary (virtual) coordinate field information

< Pointer to the ray identifier information


Now that the particles have been created, they can be iterated over and each marched one step in space. The global indices of the local ray segment storage can be easily accessed and appended. This forms a local collection of globally index ray segments.

< Recalculate the number of particles that are in the domain

< Count the number of steps that the particles have taken

< WHILE THERE ARE PARTICLES IN ANY DOMAIN

Use the ParticleStep function to calculate the path lengths of the rays through each cell so that they can be stored. This function also sets up the solve particle infrastructure.

Get all of the ray information from the particle Get the ntheta and nphi from the particle that is currently being looked at. This will be used to identify its ray and calculate its direction.

< Iterate over the particles present in the domain. How to isolate the particles in this domain and iterate over them? If there are no

IF THE CELL NUMBER IS RETURNED NEGATIVE, THEN WE HAVE REACHED THE BOUNDARY OF THE DOMAIN >> This exits the loop This function returns multiple values if multiple points are input to it Make sure that whatever cell is returned is in the radiation domain Assemble a vector of vectors etc associated with each cell index, angular coordinate, and space step? The boundary has been reached if any of these conditions don't hold

Step 3.5: The cells need to be removed if they are inside a boundary cell. Therefore, after each step (where the particle location is fed in) check for whether the cell is still within the interior region. If it is not (if it's in a boundary cell) then it should be deleted here. Condition for one dimensional domains to avoid infinite rays perpendicular to the x-axis If the domain is 1D and the x-direction of the particle is zero then delete the particle here

If the boundary has been reached by this ray, then add a boundary condition segment to the ray.

Delete the search particle associated with the ray

< Delete the particle!

< Check the point replacing the one that was deleted

Step 4: Push the particle virtual coordinates to the intersection that was found in the previous step. This ensures that the next calculated path length will start from the boundary of the adjacent cell.

< Only use the literal intersection coordinate if it exists. This will be decided above.

Step 5: Instead of using the cell face to step into the opposite cell, step the physical coordinates just beyond the intersection. This avoids issues with hitting corners and potential ghost cell weirdness. It will be slower than the face flipping but it will be more reliable. Update the coordinates of the particle. It doesn't matter which method is used, this will be the same procedure.

< Update the coordinates of the particle to move it beyond the face of the adjacent cell.

< Reset the path length to zero

Restore the fields associated with the particles after all of the particles have been stepped

DMSwarm Migrate to move the ray search particle into the next domain if it has crossed. If it no longer appears in this domain then end the ray segment.

< Migrate the search particles and remove the particles that have left the domain space.

< Update the loop condition. Recalculate the number of particles that are in the domain.

< Update the loop condition. Recalculate the number of particles that are in the domain.

< Pointer to the ray identifier information

< Get the fields from the radsolve swarm so the new point can be written to them

< Get the fields from the radsolve swarm so the new point can be written to them

Size each of the entries to hold all of the wavelengths being transported.

= 2 * (the number of independant wavelengths that are being considered). Should be read from absorption model.

Reimplemented in ablate::radiation::SurfaceRadiation.

◆ ParticleStep()

void ablate::radiation::Radiation::ParticleStep ( ablate::domain::SubDomain subDomain,
DM  faceDM,
const PetscScalar *  faceGeomArray,
DM  radReturn,
PetscInt  nlocalpoints,
PetscInt  nglobalpoints 
)
virtual

Routine to move the particle one step.

Determines the next location of the search particles during the initialization

Check that the particle is in a valid region

Declare some information associated with the field declarations

< Pointer to the primary (virtual) coordinate field information

< Pointer to the ray identifier information

Get all of the ray information from the particle Get the ntheta and nphi from the particle that is currently being looked at. This will be used to identify its ray and calculate its direction.

Check that the particle is in a valid region


The face stepping routine will give the precise path length of the mesh without any error. It will also allow the faces of the cells to be accounted for so that the boundary conditions and the conditions at reflection can be accounted for. This will make the entire initialization much faster by only requiring a single step through each cell. Additionally, the option for reflection is opened because the faces and their normals are now more easily accessed during the initialization. In the future, the carrier particles may want to be given some information that the boundary label carries when the search particle happens upon it so that imperfect reflection can be implemented.

Step 1: Register the current cell index in the rays vector. The physical coordinates that have been set in the previous step / loop will be immediately registered. Because the ray comes from the origin, all of the cell indexes are naturally ordered from the center out

Step 2: Acquire the intersection of the particle search line with the segment or face. In the case if a two dimensional mesh, the virtual coordinate in the z direction will need to be solved for because the three dimensional line will not have a literal intersection with the segment of the cell. The third coordinate can be solved for in this case. Here we are figuring out what distance the ray spends inside the cell that it has just registered.

March over each face on this cell in order to check them for the one which intersects this ray next

< Get the face geometry associated with the current cell

Check every face for intersection with the segment. The segment with the shortest path length for intersection will be the one that physically intercepts with the cell face and not with the nonphysical plane beyond the face.

< Reads the cell location from the current cell

Get the intersection of the direction vector with the cell face Use the plane equation and ray segment equation in order to get the face intersection with the shortest path length This will be the next position of the search particle

< Use plane intersection equation by getting the centroid and normal vector of the face

Step 3: Take this path if it is shorter than the previous one, getting the shortest path. The path should never be zero if the forwardIntersect check is functioning properly.

< Dumb check to ensure that the path length is always updated

Get the shortest path length of all of the faces. The point must be in the direction that the ray is travelling in order to be valid.

Reimplemented in ablate::radiation::RaySharingRadiation.

◆ Setup()

void ablate::radiation::Radiation::Setup ( const ablate::domain::Range cellRange,
ablate::domain::SubDomain subDomain 
)
virtual

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

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 in ablate::radiation::RaySharingRadiation, and ablate::radiation::OrthogonalRadiation.

◆ SurfaceComponent()

PetscReal ablate::radiation::Radiation::SurfaceComponent ( const PetscReal  normal[],
PetscInt  iCell,
PetscInt  nphi,
PetscInt  ntheta 
)
virtual

Determines what component of the incoming radiation should be accounted for when evaluating the irradiation for each ray. Dummy function that doesn't do anything unless it is overridden by the surface implementation

Reimplemented in ablate::radiation::SurfaceRadiation.

◆ UpdateCoordinates()

void ablate::radiation::Radiation::UpdateCoordinates ( PetscInt  ipart,
Virtualcoord virtualcoord,
PetscReal *  coord,
PetscReal  adv 
) const
protected

Update the coordinates of the particle using the virtual coordinates Moves the particle in physical space instead of only updating the virtual coordinates This function must be run on every updated particle before swarm migrate is used

Parameters
ipartthe particle index which is being updated
virtualcoordthe virtual coordinate field which is being read from
coordthe DMSwarm coordinate field which is being written to
adva multiple of the minimum cell radius by which to advance the DMSwarm coordinates ahead of the virtual coordinates

< If there are only two dimensions in this simulation

< Update the two physical coordinates

< If there are three dimensions in this simulation

< Update the three physical coordinates

Member Data Documentation

◆ dim

PetscInt ablate::radiation::Radiation::dim = 0
protected

Class inputs and Variables.

Number of dimensions that the domain exists within


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