ABLATE Source Documentation  0.12.34
raySharingRadiation.hpp
1 #ifndef ABLATELIBRARY_RAYSHARINGRADIATION_H
2 #define ABLATELIBRARY_RAYSHARINGRADIATION_H
3 
4 #include "domain/reverseRange.hpp"
5 #include "radiation.hpp"
6 
7 namespace ablate::radiation {
8 
10  public:
11  RaySharingRadiation(const std::string& solverId, const std::shared_ptr<domain::Region>& region, const PetscInt raynumber,
12  std::shared_ptr<eos::radiationProperties::RadiationModel> radiationModelIn, std::shared_ptr<ablate::monitors::logs::Log> = {});
14 
20  void Setup(const ablate::domain::Range& cellRange, ablate::domain::SubDomain& subDomain) override;
21 
28  void IdentifyNewRaysOnRank(ablate::domain::SubDomain& subDomain, DM radReturn, PetscInt npoints) override;
29 
40  void ParticleStep(ablate::domain::SubDomain& subDomain, DM faceDM, const PetscScalar* faceGeomArray, DM radReturn, PetscInt nlocalpoints,
41  PetscInt nglobalpoints) override;
42 
43  static inline std::string GetClassType() { return "RaySharingRadiation"; }
44 
52  void SetBoundary(CellSegment& raySegment, PetscInt index, Identifier identifier) override {
53  PetscMPIInt rank;
54  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
55  if (identifier.originRank == rank) {
56  raySegment.cell = index;
57  raySegment.pathLength = -1;
58  }
59  }
60 
61  protected:
64 
66  std::vector<PetscReal> remoteMap;
67 };
68 
69 } // namespace ablate::radiation
70 
71 #endif // ABLATELIBRARY_RAYSHARINGRADIATION_H
Definition: subDomain.hpp:19
Definition: radiation.hpp:25
std::string solverId
the name of this solver
Definition: radiation.hpp:258
const std::shared_ptr< domain::Region > region
the region for which this solver applies
Definition: radiation.hpp:261
Definition: raySharingRadiation.hpp:9
void SetBoundary(CellSegment &raySegment, PetscInt index, Identifier identifier) override
Definition: raySharingRadiation.hpp:52
void IdentifyNewRaysOnRank(ablate::domain::SubDomain &subDomain, DM radReturn, PetscInt npoints) override
Definition: raySharingRadiation.cpp:17
void Setup(const ablate::domain::Range &cellRange, ablate::domain::SubDomain &subDomain) override
Definition: raySharingRadiation.cpp:9
ablate::domain::ReverseRange indexLookup
used to look up from the cell id to range index
Definition: raySharingRadiation.hpp:63
void ParticleStep(ablate::domain::SubDomain &subDomain, DM faceDM, const PetscScalar *faceGeomArray, DM radReturn, PetscInt nlocalpoints, PetscInt nglobalpoints) override
Routine to move the particle one step.
Definition: raySharingRadiation.cpp:78
std::vector< PetscReal > remoteMap
the indexes mapping to the ray id
Definition: raySharingRadiation.hpp:66
Definition: range.hpp:11
Definition: reverseRange.hpp:11
Definition: radiation.hpp:172
PetscInt cell
< Stores the cell indices of the segment locally.
Definition: radiation.hpp:174
Definition: radiation.hpp:44
PetscInt originRank
the rank for the start of the ray
Definition: radiation.hpp:46