ABLATE Source Documentation  0.12.33
domain.hpp
1 #ifndef ABLATELIBRARY_DOMAIN_H
2 #define ABLATELIBRARY_DOMAIN_H
3 #include <petsc.h>
4 #include <algorithm>
5 #include <map>
6 #include <memory>
7 #include <string>
8 #include <utility>
9 #include <vector>
10 #include "domain/fieldDescriptor.hpp"
11 #include "domain/modifiers/modifier.hpp"
12 #include "fieldDescription.hpp"
13 #include "initializer.hpp"
14 #include "io/serializable.hpp"
15 #include "mathFunctions/fieldFunction.hpp"
16 #include "region.hpp"
17 #include "utilities/loggable.hpp"
18 #include "utilities/nonCopyable.hpp"
19 
20 namespace ablate::solver {
21 // forward declare the Solver
22 class Solver;
23 } // namespace ablate::solver
24 
25 namespace ablate::domain {
26 // forward declare the subDomain
27 class SubDomain;
28 
29 class Domain : private utilities::Loggable<Domain>, private ablate::utilities::NonCopyable {
30  public:
32  const static inline std::string solution_vector_name = "solution";
33 
35  const static inline std::string aux_vector_name = "aux";
36 
37  protected:
38  Domain(DM dm, std::string name, std::vector<std::shared_ptr<FieldDescriptor>>, std::vector<std::shared_ptr<modifiers::Modifier>> modifiers,
39  const std::shared_ptr<parameters::Parameters>& options = {}, bool setFromOptions = true);
40  virtual ~Domain();
41 
42  // The primary dm
43  DM dm;
44 
45  private:
46  // the name of the dm
47  std::string name;
48 
49  // Hold a copy of the comm for this DM
50  const MPI_Comm comm;
51 
52  // List of classes that are used to describe fields
53  const std::vector<std::shared_ptr<FieldDescriptor>> fieldDescriptors;
54 
55  // Keep track of all solution fields
56  std::vector<Field> fields;
57 
58  // This domain can be partitions into multiple subdomains
59  std::vector<std::shared_ptr<SubDomain>> subDomains;
60 
62  Vec solGlobalField;
63 
64  void CreateStructures();
65 
66  // keep a list of functions that modify the dm
67  std::vector<std::shared_ptr<modifiers::Modifier>> modifiers;
68 
70  PetscOptions petscOptions = nullptr;
71 
72  public:
73  [[nodiscard]] const std::string& GetName() const { return name; }
74 
75  inline DM& GetDM() noexcept { return dm; }
76 
81  inline Vec GetSolutionVector() { return solGlobalField; }
82 
87  void RegisterField(const ablate::domain::FieldDescription& fieldDescription);
88 
89  [[nodiscard]] PetscInt GetDimensions() const noexcept;
90 
96  void InitializeSubDomains(const std::vector<std::shared_ptr<solver::Solver>>& solvers = {}, const std::shared_ptr<ablate::domain::Initializer>& initializations = {},
97  const std::vector<std::shared_ptr<mathFunctions::FieldFunction>>& = {});
98 
104  void ProjectFieldFunctions(const std::vector<std::shared_ptr<mathFunctions::FieldFunction>>& fieldFunctions, Vec globVec, PetscReal time = 0.0);
105 
111  std::shared_ptr<SubDomain> GetSubDomain(const std::shared_ptr<Region>& name) const;
112 
116  std::vector<std::weak_ptr<io::Serializable>> GetSerializableSubDomains();
117 
123  [[nodiscard]] inline const Field& GetField(int fieldId) const { return fields[fieldId]; }
124 
125  [[nodiscard]] inline const Field& GetField(const std::string_view& fieldName) const {
126  auto field = std::find_if(fields.begin(), fields.end(), [&fieldName](auto field) { return field.name == fieldName; });
127  if (field != fields.end()) {
128  return *field;
129  } else {
130  throw std::invalid_argument("Cannot locate field with name " + std::string(fieldName) + " in domain " + name);
131  }
132  }
133 
139  [[nodiscard]] inline const std::vector<Field>& GetFields() const { return fields; }
140 
146  bool CheckFieldValues(Vec globSourceVector = nullptr);
147 };
148 } // namespace ablate::domain
149 #endif // ABLATELIBRARY_DOMAIN_H
Definition: domain.hpp:29
bool CheckFieldValues(Vec globSourceVector=nullptr)
Definition: domain.cpp:252
void InitializeSubDomains(const std::vector< std::shared_ptr< solver::Solver >> &solvers={}, const std::shared_ptr< ablate::domain::Initializer > &initializations={}, const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> &={})
Definition: domain.cpp:143
void ProjectFieldFunctions(const std::vector< std::shared_ptr< mathFunctions::FieldFunction >> &fieldFunctions, Vec globVec, PetscReal time=0.0)
Definition: domain.cpp:198
void RegisterField(const ablate::domain::FieldDescription &fieldDescription)
Definition: domain.cpp:71
std::vector< std::weak_ptr< io::Serializable > > GetSerializableSubDomains()
Definition: domain.cpp:190
const Field & GetField(int fieldId) const
Definition: domain.hpp:123
static const std::string solution_vector_name
The name of the solution field vector.
Definition: domain.hpp:32
static const std::string aux_vector_name
The name of the aux field vector.
Definition: domain.hpp:35
const std::vector< Field > & GetFields() const
Definition: domain.hpp:139
std::shared_ptr< SubDomain > GetSubDomain(const std::shared_ptr< Region > &name) const
Definition: domain.cpp:123
Vec GetSolutionVector()
Definition: domain.hpp:81
Definition: loggable.hpp:9
Definition: nonCopyable.hpp:9
Definition: fieldDescription.hpp:20
Definition: field.hpp:17