ABLATE Source Documentation  0.12.34
swarmAccessor.hpp
1 #ifndef ABLATELIBRARY_SWARMACCESSOR_HPP
2 #define ABLATELIBRARY_SWARMACCESSOR_HPP
3 
4 #include <petsc.h>
5 #include <map>
6 #include "accessor.hpp"
7 #include "particles/field.hpp"
8 #include "utilities/petscUtilities.hpp"
9 
10 namespace ablate::particles::accessors {
14 class SwarmAccessor : public Accessor<const PetscReal> {
15  private:
17  const DM& swarmDm;
18 
20  const std::map<std::string, Field>& fieldsMap;
21 
23  Vec solutionVec;
24 
26  const PetscScalar* solutionValues{};
27 
28  public:
29  SwarmAccessor(bool cachePointData, const DM& swarmDm, const std::map<std::string, Field>& fieldsMap, Vec solutionVec)
30  : Accessor(cachePointData), swarmDm(swarmDm), fieldsMap(fieldsMap), solutionVec(solutionVec) {
31  // extract the array from the vector
32  VecGetArrayRead(solutionVec, &solutionValues) >> utilities::PetscUtilities::checkError;
33  }
34 
35  ~SwarmAccessor() override { VecRestoreArrayRead(solutionVec, &solutionValues) >> utilities::PetscUtilities::checkError; }
36 
41  [[nodiscard]] inline PetscInt GetNumberParticles() const {
42  PetscInt size;
43  DMSwarmGetLocalSize(swarmDm, &size) >> utilities::PetscUtilities::checkError;
44  return size;
45  }
46 
50  SwarmAccessor(const SwarmAccessor&) = delete;
51 
52  protected:
58  ConstPointData CreateData(const std::string& fieldName) override {
59  const auto& field = fieldsMap.at(fieldName);
60  if (field.location == domain::FieldLocation::SOL) {
61  return {solutionValues, field};
62  } else {
63  // get the field from the dm
64  PetscScalar* values;
65  DMSwarmGetField(swarmDm, field.name.c_str(), nullptr, nullptr, (void**)&values) >> utilities::PetscUtilities::checkError;
66 
67  // Register the cleanup
69  const std::string name = field.name;
70  DMSwarmRestoreField(swarmDm, name.c_str(), nullptr, nullptr, (void**)&values) >> utilities::PetscUtilities::checkError;
71  });
72 
73  return {values, field};
74  }
75  }
76 };
77 } // namespace ablate::particles::accessors
78 #endif // ABLATELIBRARY_SWARMACCESSOR_HPP
Definition: accessor.hpp:17
void RegisterCleanupFunction(const std::function< void()> &function)
Definition: accessor.hpp:78
Definition: swarmAccessor.hpp:14
SwarmAccessor(const SwarmAccessor &)=delete
PetscInt GetNumberParticles() const
Definition: swarmAccessor.hpp:41
ConstPointData CreateData(const std::string &fieldName) override
Definition: swarmAccessor.hpp:58
static utilities::PetscUtilities::ErrorChecker checkError
Definition: petscUtilities.hpp:46
Definition: pointData.hpp:14