Chi-Tech
ff_gridbased_00_constrdestr.cc
Go to the documentation of this file.
2
8
11
12#include "ChiObjectFactory.h"
13
14namespace chi_physics
15{
16
18
20{
22
23 params.SetDocGroup("DocFieldFunction");
24
25 params.AddRequiredParameter<std::string>(
26 "sdm_type", "The spatial discretization type to be used");
27
29 "initial_value", 0.0, "The initial value to assign to the field function");
30
32 "quadrature_order",
33 0,
34 "If supplied, will overwrite the default for the "
35 "specific discretization-coordinate system combination.");
36
37 params.AddOptionalParameter("coordinate_system",
38 "cartesian",
39 "Coordinate system to apply to element mappings");
40
41 using namespace chi_data_types;
43 "sdm_type",
44 AllowableRangeList::New({"FV", "PWLC", "PWLD", "LagrangeC", "LagrangeD"}));
46 "coordinate_system",
47 AllowableRangeList::New({"cartesian", "cylindrical", "spherical"}));
48
49 return params;
50}
51
52// ##################################################################
53/**ObjectMaker based constructor.*/
55 const chi::InputParameters& params)
56 : FieldFunction(params),
57 sdm_(MakeSpatialDiscretization(params)),
58 ghosted_field_vector_(MakeFieldVector(*sdm_, GetUnknownManager())),
59 local_grid_bounding_box_(
60 chi_mesh::GetCurrentHandler().GetGrid()->GetLocalBoundingBox())
61{
62 ghosted_field_vector_->Set(params.GetParamValue<double>("initial_value"));
63}
64
65// ##################################################################
67 const std::string& text_name,
68 chi_math::SDMPtr& discretization_ptr,
69 chi_math::Unknown unknown)
70 : FieldFunction(text_name, std::move(unknown)),
71 sdm_(discretization_ptr),
72 ghosted_field_vector_(MakeFieldVector(*sdm_, GetUnknownManager())),
73 local_grid_bounding_box_(sdm_->Grid().GetLocalBoundingBox())
74{
75}
76
77// ##################################################################
79 const std::string& text_name,
80 chi_math::SDMPtr& sdm_ptr,
81 chi_math::Unknown unknown,
82 const std::vector<double>& field_vector)
83 : FieldFunction(text_name, std::move(unknown)),
84 sdm_(sdm_ptr),
85 ghosted_field_vector_(MakeFieldVector(*sdm_, GetUnknownManager())),
86 local_grid_bounding_box_(sdm_->Grid().GetLocalBoundingBox())
87{
89 field_vector.size() != ghosted_field_vector_->LocalSize(),
90 "Constructor called with incompatible size field vector.");
91
92 ghosted_field_vector_->Set(field_vector);
93}
94
95// ##################################################################
97 chi_math::SDMPtr& sdm_ptr,
98 chi_math::Unknown unknown,
99 double field_value)
100 : FieldFunction(text_name, std::move(unknown)),
101 sdm_(sdm_ptr),
102 ghosted_field_vector_(MakeFieldVector(*sdm_, GetUnknownManager())),
103 local_grid_bounding_box_(sdm_->Grid().GetLocalBoundingBox())
104{
105 ghosted_field_vector_->Set(field_value);
106}
107
108// ##################################################################
109/**Returns the spatial discretization method.*/
112{
113 return *sdm_;
114}
115
116// ##################################################################
117/**Returns a read-only reference to the locally stored field data.*/
118const std::vector<double>& FieldFunctionGridBased::FieldVectorRead() const
119{
120 return ghosted_field_vector_->LocalSTLData();
121}
122/**Returns a reference to the locally stored field data.*/
124{
125 return ghosted_field_vector_->LocalSTLData();
126}
127
128// ##################################################################
129/**Private method for creating the spatial discretization method.*/
131 const chi::InputParameters& params)
132{
133 const auto& user_params = params.ParametersAtAssignment();
134 const auto& grid_ptr = chi_mesh::GetCurrentHandler().GetGrid();
135 const auto sdm_type = params.GetParamValue<std::string>("sdm_type");
136
142
143 if (sdm_type == "FV") return FV::New(*grid_ptr);
144
147 std::string cs = "cartesian";
148 if (user_params.Has("coordinate_system"))
149 {
150 cs = params.GetParamValue<std::string>("coordinate_system");
151
152 using namespace chi_math;
153 if (cs == "cartesian") cs_type = CoordinateSystemType::CARTESIAN;
154 if (cs == "cylindrical") cs_type = CoordinateSystemType::CYLINDRICAL;
155 if (cs == "spherical") cs_type = CoordinateSystemType::SPHERICAL;
156 }
157
159
160 if (user_params.Has("quadrature_order"))
161 {
162 const uint32_t max_order =
163 static_cast<uint32_t>(chi_math::QuadratureOrder::FORTYTHIRD);
164 const uint32_t q_order_int =
165 params.GetParamValue<uint32_t>("quadrature_order");
166 ChiInvalidArgumentIf(q_order_int > max_order,
167 "Invalid quadrature point order " +
168 std::to_string(q_order_int));
169 q_order = static_cast<chi_math::QuadratureOrder>(q_order_int);
170 }
171 else // Defaulted
172 {
173 if (cs == "cartesian") q_order = chi_math::QuadratureOrder::SECOND;
174 if (cs == "cylindrical") q_order = chi_math::QuadratureOrder::THIRD;
175 if (cs == "spherical") q_order = chi_math::QuadratureOrder::FOURTH;
176 }
177
178 if (sdm_type == "PWLC") return PWLC::New(*grid_ptr, q_order, cs_type);
179 else if (sdm_type == "PWLD")
180 return PWLD::New(*grid_ptr, q_order, cs_type);
181 else if (sdm_type == "LagrangeC")
182 return LagC::New(*grid_ptr, q_order, cs_type);
183 else if (sdm_type == "LagrangeD")
184 return LagD::New(*grid_ptr, q_order, cs_type);
185
186 // If not returned by now
187 ChiInvalidArgument("Unsupported sdm_type \"" + sdm_type + "\"");
188}
189
190/**Private method for creating the field vector.*/
191std::unique_ptr<chi_math::GhostedParallelSTLVector>
193 const chi_math::SpatialDiscretization& discretization,
194 const chi_math::UnknownManager& uk_man)
195{
196 auto field = std::make_unique<chi_math::GhostedParallelSTLVector>(
197 discretization.GetNumLocalDOFs(uk_man),
198 discretization.GetNumGlobalDOFs(uk_man),
199 discretization.GetGhostDOFIndices(uk_man),
200 Chi::mpi.comm);
201
202 return field;
203}
204
205} // namespace chi_physics
#define ChiInvalidArgumentIf(condition, message)
#define ChiInvalidArgument(message)
static chi::MPI_Info & mpi
Definition: chi_runtime.h:78
void SetDocGroup(const std::string &doc_group)
void AddRequiredParameter(const std::string &name, const std::string &doc_string)
const ParameterBlock & ParametersAtAssignment() const
void ConstrainParameterRange(const std::string &param_name, AllowableRangePtr allowable_range)
void AddOptionalParameter(const std::string &name, T value, const std::string &doc_string)
const MPI_Comm & comm
MPI communicator.
Definition: mpi_info.h:28
T GetParamValue(const std::string &param_name) const
virtual std::vector< int64_t > GetGhostDOFIndices(const UnknownManager &unknown_manager) const =0
size_t GetNumGlobalDOFs(const UnknownManager &unknown_manager) const
size_t GetNumLocalDOFs(const UnknownManager &unknown_manager) const
chi_mesh::MeshContinuumPtr & GetGrid() const
static std::unique_ptr< chi_math::GhostedParallelSTLVector > MakeFieldVector(const chi_math::SpatialDiscretization &discretization, const chi_math::UnknownManager &uk_man)
const std::vector< double > & FieldVectorRead() const
FieldFunctionGridBased(const chi::InputParameters &params)
static chi::InputParameters GetInputParameters()
const chi_math::SpatialDiscretization & GetSpatialDiscretization() const
std::unique_ptr< chi_math::GhostedParallelSTLVector > ghosted_field_vector_
static chi_math::SDMPtr MakeSpatialDiscretization(const chi::InputParameters &params)
static chi::InputParameters GetInputParameters()
QuadratureOrder
Definition: quadrature.h:12
std::shared_ptr< SpatialDiscretization > SDMPtr
CoordinateSystemType
Definition: chi_math.h:29
MeshHandler & GetCurrentHandler()
RegisterChiObject(chi_physics, FieldFunctionGridBased)