v0.15.0
Loading...
Searching...
No Matches
FieldEvaluator.hpp
Go to the documentation of this file.
1/** \file FieldEvaluator.hpp
2 * \brief Field Evaluator
3 *
4 * Evaluate field at given coordinate
5 *
6 *
7 * \ingroup field_evaluator
8 */
9
10#ifndef __FIELD_EVALUATOR_HPP__
11 #define __FIELD_EVALUATOR_HPP__
12
13 #include "UnknownInterface.hpp"
14
15namespace MoFEM {
16
17/** \brief Field evaluator interface
18
19 * \ingroup field_evaluator
20 */
22
25
26 /**
27 * @brief Default evaluator for setting integration points
28 *
29 */
31
32 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
33 UnknownInterface **iface) const;
34
36
37 /**
38 * @brief Get the Data object
39 *
40 * Pack pointers with data structures for field evaluator and finite
41 * element. Function return shared pointer if returned shared pointer
42 * is reset; all data are destroyed. The idea is to pack all data in
43 * one structure, create shared pointer to it and return aliased shared
44 * pointer to one of the elements of the data structure. It is a bit
45 * complicated, but has great flexibility.
46 *
47 * @tparam VE
48 * @tparam SetIntegrationPtsMethodData
49 * @tparam SetIntegrationPtsMethod
50 * @param ptr
51 * @param nb_eval_points
52 * @param eps
53 * @param verb
54 * @return boost::shared_ptr<SPD>
55 */
56 template <typename VE, typename SPD = SetIntegrationPtsMethodData,
57 typename SP = SetIntegrationPtsMethod>
58 boost::shared_ptr<SPD>
59 getData(const double *ptr = nullptr, const int nb_eval_points = 0,
60 const double eps = 1e-12, VERBOSITY_LEVELS verb = QUIET);
61
62 /**
63 * @brief Build spatial tree
64 *
65 * @param finite_element finite element name
66 * @return MoFEMErrorCode
67 */
68 template <int D>
69 inline MoFEMErrorCode
70 buildTree(boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
71 const std::string finite_element) {
72 return buildTree<D>(spd_ptr, finite_element, BitRefLevel(), BitRefLevel());
73 }
74
75 /**
76 * @brief Build spatial tree
77 *
78 * @param finite_element finite element name
79 * @param bit bit reference level
80 * @param mask mask reference level
81 * @return MoFEMErrorCode
82 */
83 template <int D>
85 buildTree(boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
86 const std::string finite_element, BitRefLevel bit,
87 BitRefLevel mask);
88
89 /**
90 * @brief Evaluate field at actbitray position
91 *
92 * \code
93
94 std::array<double, 3> point = {0, 0, 0};
95 const double dist = 0.3;
96 std::array<double, 6> eval_points = {-1., -1., -1., 1., 1., 1. };
97
98 using VolEle = VolumeElementForcesAndSourcesCore;
99 auto data_ptr =
100 m_field.getInterface<FieldEvaluatorInterface>()->
101 getData<VolEle>(point.data(), point.size/3);
102
103 if(auto vol_ele = data_ptr->feMethod.lock()) {
104 // push operators to finite element instance, e.g.
105 vol_ele->getOpPtrVector().push_back(new MyOp());
106 // iterate over elemnts with evaluated points
107 auto cache_ptr = boost::make_shared<CacheTuple>();
108 CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
109 // SPACE_DIM is 3 in this example
110 CHKERR m_field.getInterface<FieldEvaluatorInterface>()
111 ->evalFEAtThePoint<3>(point.data(), dist, prb_ptr->getName(),
112 "FINITE_ELEMENT_NAME",
113 data_ptr, m_field.get_comm_rank(),
114 m_field.get_comm_rank(), cache_ptr);
115 }
116
117 * \endcode
118 *
119 * @param point point used to find tetrahedrons
120 * @param distance distance from the point where tetrahedrons are searched
121 * @param problem problem name
122 * @param finite_element finite element name
123 * @param data_ptr pointer to data abut gauss points
124 * @param lower_rank lower processor
125 * @param upper_rank upper process
126 * @param cache_ptr cache or problem entities
127 * @param bh control looping over entities, e.g. throwing error if element not
128 found
129 * @param verb verbosity level
130 * @return MoFEMErrorCode
131 */
132 template <int D>
133 inline MoFEMErrorCode
134 evalFEAtThePoint(const double *const point, const double distance,
135 const std::string problem, const std::string finite_element,
136 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
137 int lower_rank, int upper_rank,
138 boost::shared_ptr<CacheTuple> cache_ptr,
140 return evalFEAtThePoint<D>(point, distance, problem, finite_element,
141 data_ptr, lower_rank, upper_rank, BitRefLevel(),
142 BitRefLevel(), cache_ptr, bh, verb);
143 }
144
145 using VecDataPts = std::vector<
146
147 std::pair<boost::shared_ptr<VectorDouble>,
148 boost::shared_ptr<VectorDouble>>
149
150 >;
151
152 using MatDataPts = std::vector<
153
154 std::pair<boost::shared_ptr<MatrixDouble>,
155 boost::shared_ptr<MatrixDouble>>
156
157 >;
158
159 /**
160 * @brief Evaluate field values at integration points of finite elements
161 * to which operator is pushed
162 *
163 * @param vec_data_pts vector of data shared between physical and post
164 * proc elements
165 * @param mat_data_pts matrix of data shared between physical and post
166 * proc elements
167 * @param finite_element finite element name on which field is evaluated
168 * @param data_ptr data about field evaluator
169 * @param bh control looping over entities, e.g. throwing error if element not
170 * found
171 *
172 */
173 template <int D>
175 getDataOperator(VecDataPts vec_data_pts, MatDataPts mat_data_pts,
176 const std::string finite_element,
177 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
178 int lower_rank, int upper_rank, BitRefLevel bit,
180 VERBOSITY_LEVELS verb = QUIET);
181
182 /** @brief Evaluate field at actbitray position
183 *
184 * @param point point used to find tetrahedrons
185 * @param distance distance from the point where tetrahedrons are searched
186 * @param problem problem name
187 * @param finite_element finite element name
188 * @param data_ptr pointer to data abut gauss points
189 * @param lower_rank lower processor
190 * @param upper_rank upper process
191 * @param bit bit reference level
192 * @param mask mask reference level
193 * @param cache_ptr cache or problem entities
194 * @param bh control looping over entities, e.g. throwing error if element not
195 found
196 * @param verb verbosity level
197 * @return MoFEMErrorCode
198 */
199 template <int D>
201 evalFEAtThePoint(const double *const point, const double distance,
202 const std::string problem, const std::string finite_element,
203 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
204 int lower_rank, int upper_rank, BitRefLevel bit,
205 BitRefLevel mask, boost::shared_ptr<CacheTuple> cache_ptr,
207
208private:
210
211 template <int D>
213 buildTreeImpl(boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
214 const std::string finite_element, BitRefLevel bit,
215 BitRefLevel mask);
216
217 template <int D>
219 const double *const point, const double distance,
220 const std::string problem, const std::string finite_element,
221 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
222 int upper_rank, BitRefLevel bit, BitRefLevel mask,
223 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh = MF_EXIST,
224 VERBOSITY_LEVELS verb = QUIET);
225
226 template <int D>
228 getDataOperatorImpl(VecDataPts vec_data_pts, MatDataPts mat_data_pts,
229 const std::string finite_element,
230 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
231 int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask,
233};
234
235template <>
236MoFEMErrorCode FieldEvaluatorInterface::buildTree<2>(
237 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
238 const std::string finite_element, BitRefLevel bit, BitRefLevel mask);
239
240template <>
241MoFEMErrorCode FieldEvaluatorInterface::buildTree<3>(
242 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
243 const std::string finite_element, BitRefLevel bit, BitRefLevel mask);
244
245template <>
246MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<2>(
247 const double *const point, const double distance, const std::string problem,
248 const std::string finite_element,
249 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
250 int upper_rank, BitRefLevel bit, BitRefLevel mask,
251 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
252 VERBOSITY_LEVELS verb);
253
254template <>
255MoFEMErrorCode FieldEvaluatorInterface::evalFEAtThePoint<3>(
256 const double *const point, const double distance, const std::string problem,
257 const std::string finite_element,
258 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
259 int upper_rank, BitRefLevel bit, BitRefLevel mask,
260 boost::shared_ptr<CacheTuple> cache_ptr, MoFEMTypes bh,
261 VERBOSITY_LEVELS verb);
262
263template <>
264ForcesAndSourcesCore::UserDataOperator *FieldEvaluatorInterface::getDataOperator<2>(
267 const std::string finite_element,
268 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
269 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
270 VERBOSITY_LEVELS verb);
271
272template <>
273ForcesAndSourcesCore::UserDataOperator *FieldEvaluatorInterface::getDataOperator<3>(
276 const std::string finite_element,
277 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
278 int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh,
279 VERBOSITY_LEVELS verb);
280
281
283
285
286 /**
287 * @brief Set the Gauss Pts data
288 *
289 * @param fe_method_ptr pointer to finite element instance
290 * @param eval_points pointer to array with evaluation points
291 * @param nb_eval_points number of evaluated points
292 * @param eps tolerance used to find if point is in the element
293 * @param verb
294 */
296 boost::shared_ptr<MoFEM::ForcesAndSourcesCore> fe_method_ptr,
297 const double *eval_points, const int nb_eval_points, const double eps,
299
300 void setEvalPoints(const double *ptr, const int nb_eval_points);
301
302 boost::weak_ptr<MoFEM::ForcesAndSourcesCore> feMethodPtr;
304 boost::scoped_ptr<AdaptiveKDTree> treePtr;
305
306protected:
307 const double *evalPoints;
309 double eps;
311
314 std::vector<EntityHandle> evalPointEntityHandle;
315
318};
319
323 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr);
324 MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row,
325 int order_col, int order_data);
326
327private:
328 boost::weak_ptr<SetIntegrationPtsMethodData> dataPtr;
329};
330
331template <typename VE, typename SPD, typename SP>
332boost::shared_ptr<SPD>
333FieldEvaluatorInterface::getData(const double *ptr, const int nb_eval_points,
334 const double eps, VERBOSITY_LEVELS verb) {
335
336 struct PackData {
337 boost::scoped_ptr<VE> elePtr; //< finite element pointer
338 boost::scoped_ptr<SPD> setPtsDataPtr; //< integration pts data
339 boost::scoped_ptr<SP> setPtsPtr; //< integration pts method
340 };
341 boost::shared_ptr<PackData> pack_data(new PackData());
342 MoFEM::Interface &m_field = cOre;
343 pack_data->elePtr.reset(new VE(m_field));
344
345 pack_data->setPtsDataPtr.reset(
346 new SPD(boost::shared_ptr<VE>(pack_data, pack_data->elePtr.get()), ptr,
347 nb_eval_points, eps, verb));
348 pack_data->setPtsPtr.reset(new SP(
349 boost::shared_ptr<SPD>(pack_data, pack_data->setPtsDataPtr.get())));
350 pack_data->elePtr->setRuleHook = boost::ref(*pack_data->setPtsPtr);
351 boost::shared_ptr<SPD> data(pack_data, pack_data->setPtsDataPtr.get());
352 return data;
353}
354
355} // namespace MoFEM
356
357#endif // __FIELD_EVALUATOR_HPP__
358
359/**
360 * \defgroup field_evaluator Field Evaluator
361 * \brief Evaluate field at the point
362 *
363 * \ingroup mofem
364 */
MoFEM interface.
static const double eps
VERBOSITY_LEVELS
Verbosity levels.
@ QUIET
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ MF_EXIST
auto bit
set bit
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
Core (interface) class.
Definition Core.hpp:82
Deprecated interface functions.
void setEvalPoints(const double *ptr, const int nb_eval_points)
boost::weak_ptr< MoFEM::ForcesAndSourcesCore > feMethodPtr
MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
boost::weak_ptr< SetIntegrationPtsMethodData > dataPtr
Field evaluator interface.
MoFEMErrorCode evalFEAtThePoint(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field at actbitray position.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode evalFEAtThePoint(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field at actbitray position.
ForcesAndSourcesCore::UserDataOperator * getDataOperatorImpl(VecDataPts vec_data_pts, MatDataPts mat_data_pts, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
std::vector< std::pair< boost::shared_ptr< MatrixDouble >, boost::shared_ptr< MatrixDouble > > > MatDataPts
boost::shared_ptr< SPD > getData(const double *ptr=nullptr, const int nb_eval_points=0, const double eps=1e-12, VERBOSITY_LEVELS verb=QUIET)
Get the Data object.
MoFEMErrorCode buildTreeImpl(boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
MoFEMErrorCode buildTree(boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode buildTree(boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
Build spatial tree.
ForcesAndSourcesCore::UserDataOperator * getDataOperator(VecDataPts vec_data_pts, MatDataPts mat_data_pts, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field values at integration points of finite elements to which operator is pushed.
MoFEMErrorCode evalFEAtThePointImpl(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetIntegrationPtsMethodData > data_ptr, int lower_rank, int upper_rank, BitRefLevel bit, BitRefLevel mask, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
std::vector< std::pair< boost::shared_ptr< VectorDouble >, boost::shared_ptr< VectorDouble > > > VecDataPts
structure to get information from mofem into EntitiesFieldData
base class for all interface classes