v0.13.2
Loading...
Searching...
No Matches
FieldEvaluator.cpp
Go to the documentation of this file.
1/** \file FieldEvaluator.cpp
2 * \brief Field evaluator
3 *
4 */
5
6namespace MoFEM {
7
9 : cOre(const_cast<MoFEM::Core &>(core)) {
10
11 if (!LogManager::checkIfChannelExist("FieldEvaluatorWorld")) {
12 auto core_log = logging::core::get();
13
15 "FieldEvaluatorWorld"));
17 "FieldEvaluatorSync"));
19 "FieldEvaluatorSelf"));
20
21 LogManager::setLog("FieldEvaluatorWorld");
22 LogManager::setLog("FieldEvaluatorSync");
23 LogManager::setLog("FieldEvaluatorSelf");
24
25 MOFEM_LOG_TAG("FieldEvaluatorWorld", "FieldEvaluator");
26 MOFEM_LOG_TAG("FieldEvaluatorSync", "FieldEvaluator");
27 MOFEM_LOG_TAG("FieldEvaluatorSelf", "FieldEvaluator");
28 }
29
30 MOFEM_LOG("FieldEvaluatorWorld", Sev::noisy) << "Field evaluator intreface";
31}
32
34 boost::typeindex::type_index type_index, UnknownInterface **iface) const {
36 *iface = const_cast<FieldEvaluatorInterface *>(this);
38}
39
40template <int D>
42FieldEvaluatorInterface::buildTree(boost::shared_ptr<SetPtsData> spd_ptr,
43 const std::string finite_element) {
44 CoreInterface &m_field = cOre;
46 EntityHandle fe_meshset;
47 fe_meshset = m_field.get_finite_element_meshset(finite_element);
48 Range entities;
49 CHKERR m_field.get_moab().get_entities_by_dimension(fe_meshset, D, entities,
50 true);
51 spd_ptr->treePtr.reset(new AdaptiveKDTree(&m_field.get_moab()));
52 CHKERR spd_ptr->treePtr->build_tree(entities, &spd_ptr->rooTreeSet);
54}
55
57FieldEvaluatorInterface::buildTree3D(boost::shared_ptr<SetPtsData> spd_ptr,
58 const std::string finite_element) {
59 return buildTree<3>(spd_ptr, finite_element);
60}
61
63FieldEvaluatorInterface::buildTree2D(boost::shared_ptr<SetPtsData> spd_ptr,
64 const std::string finite_element) {
65 return buildTree<2>(spd_ptr, finite_element);
66}
67
70 int order_row, int order_col,
71 int order_data) {
73
74 if (auto data_ptr = dataPtr.lock()) {
75
76 const auto nb_eval_points = data_ptr->nbEvalPoints;
77 const auto &shape_functions = data_ptr->shapeFunctions;
78 const auto &eval_pointentity_handle = data_ptr->evalPointEntityHandle;
79
80 if (auto fe_ptr = data_ptr->feMethodPtr.lock()) {
81
82 auto &fe = static_cast<ForcesAndSourcesCore &>(*fe_ptr);
83
84#ifndef NDEBUG
85 if (fe_ptr.get() != fe_raw_ptr)
86 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong FE ptr");
87#endif
88
89 const auto fe_ent = fe.numeredEntFiniteElementPtr->getEnt();
90 const auto fe_type = type_from_handle(fe_ent);
91 const auto fe_dim = moab::CN::Dimension(fe_type);
92
93 auto &gauss_pts = fe.gaussPts;
94 int nb_gauss_pts;
95
96 if (fe_dim == 3) {
97 gauss_pts.resize(4, nb_eval_points, false);
99 &shape_functions(0, 0), &shape_functions(0, 1),
100 &shape_functions(0, 2), &shape_functions(0, 3)};
102 &gauss_pts(0, 0), &gauss_pts(1, 0), &gauss_pts(2, 0)};
103
104 nb_gauss_pts = 0;
105 for (int nn = 0; nn != nb_eval_points; ++nn) {
106 if (eval_pointentity_handle[nn] == fe_ent) {
107 for (const int i : {0, 1, 2}) {
108 t_gauss_pts(i) = t_shape(i + 1);
109 }
110 gauss_pts(3, nb_gauss_pts) = nn;
111 ++t_gauss_pts;
112 ++nb_gauss_pts;
113 }
114 ++t_shape;
115 }
116 gauss_pts.resize(4, nb_gauss_pts, true);
117 } else if (fe_dim == 2) {
118 gauss_pts.resize(3, nb_eval_points, false);
120 &shape_functions(0, 0), &shape_functions(0, 1),
121 &shape_functions(0, 2)};
123 &gauss_pts(0, 0), &gauss_pts(1, 0)};
124 nb_gauss_pts = 0;
125 for (int nn = 0; nn != nb_eval_points; ++nn) {
126 if (eval_pointentity_handle[nn] == fe_ent) {
127 for (const int i : {0, 1}) {
128 t_gauss_pts(i) = t_shape(i + 1);
129 }
130 gauss_pts(2, nb_gauss_pts) = nn;
131 ++t_gauss_pts;
132 ++nb_gauss_pts;
133 }
134 ++t_shape;
135 }
136 gauss_pts.resize(3, nb_gauss_pts, true);
137 } else {
138 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
139 "Dimension not implemented");
140 }
141
142#ifndef NDEBUG
143 MOFEM_LOG("SELF", Sev::noisy)
144 << "nbEvalOPoints / nbGau_sspt_s: " << nb_eval_points << " / "
145 << nb_gauss_pts;
146 MOFEM_LOG("SELF", Sev::noisy) << "gauss pts: " << gauss_pts;
147#endif
148
149 } else
150 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
151 "Pointer to element does not exists");
152
153 } else
154 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
155 "Pointer to data does not exists");
156
158}
159
160template <int D>
162 const double *const point, const double distance, const std::string problem,
163 const std::string finite_element, boost::shared_ptr<SetPtsData> data_ptr,
164 int lower_rank, int upper_rank, boost::shared_ptr<CacheTuple> cache_ptr,
165 MoFEMTypes bh, VERBOSITY_LEVELS verb) {
166 CoreInterface &m_field = cOre;
168
169 std::vector<EntityHandle> leafs_out;
170 CHKERR data_ptr->treePtr->distance_search(point, distance, leafs_out);
171 Range tree_ents;
172 for (auto lit : leafs_out)
173 CHKERR m_field.get_moab().get_entities_by_dimension(lit, D, tree_ents,
174 true);
175 if (verb >= NOISY)
176 MOFEM_LOG("SELF", Sev::noisy) << "tree entities: " << tree_ents;
177
178 data_ptr->evalPointEntityHandle.resize(data_ptr->nbEvalPoints);
179 std::fill(data_ptr->evalPointEntityHandle.begin(),
180 data_ptr->evalPointEntityHandle.end(), 0);
181
182 std::vector<double> local_coords;
183 std::vector<double> shape;
184
185 auto nb_eval_points = data_ptr->nbEvalPoints;
186
187 for (auto tet : tree_ents) {
188 const EntityHandle *conn;
189 int num_nodes;
190 CHKERR m_field.get_moab().get_connectivity(tet, conn, num_nodes, true);
191
192 if constexpr (D == 3) {
193
194 local_coords.resize(3 * nb_eval_points);
195 shape.resize(4 * nb_eval_points);
196
197 std::array<double, 12> coords;
198 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
199
201 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
202 &local_coords[0]);
203 CHKERR Tools::shapeFunMBTET<3>(&shape[0], &local_coords[0],
204 &local_coords[1], &local_coords[2],
205 nb_eval_points);
206
209 &shape[0], &shape[1], &shape[2], &shape[3]};
211 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
212 &data_ptr->shapeFunctions(0, 2), &data_ptr->shapeFunctions(0, 3)};
213
216 &local_coords[0], &local_coords[1], &local_coords[2]};
218 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1)),
219 &(data_ptr->localCoords(0, 2))};
220
221 for (int n = 0; n != nb_eval_points; ++n) {
222
223 const double eps = data_ptr->eps;
224 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
225
226 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
227
228 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps &&
229
230 t_shape(3) >= 0 - eps && t_shape(3) <= 1 + eps) {
231
232 data_ptr->evalPointEntityHandle[n] = tet;
233 t_shape_data(i4) = t_shape(i4);
234 t_local_data(j3) = t_local(j3);
235 }
236
237 ++t_shape;
238 ++t_shape_data;
239 ++t_local;
240 ++t_local_data;
241 }
242 }
243
244 if constexpr (D == 2) {
245
246 local_coords.resize(2 * nb_eval_points);
247 shape.resize(3 * nb_eval_points);
248
249 std::array<double, 9> coords;
250 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
251
253 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
254 &local_coords[0]);
255 CHKERR Tools::shapeFunMBTRI<2>(&shape[0], &local_coords[0],
256 &local_coords[1], nb_eval_points);
257
260 &shape[0], &shape[1], &shape[2]};
262 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
263 &data_ptr->shapeFunctions(0, 2)};
264
267 &local_coords[0], &local_coords[1]};
268 data_ptr->localCoords.resize(nb_eval_points, 2);
270 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1))};
271
272 for (int n = 0; n != nb_eval_points; ++n) {
273
274 const double eps = data_ptr->eps;
275 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
276
277 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
278
279 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps) {
280
281 data_ptr->evalPointEntityHandle[n] = tet;
282 t_shape_data(i3) = t_shape(i3);
283 t_local_data(j2) = t_local(j2);
284 }
285
286 ++t_shape;
287 ++t_shape_data;
288 ++t_local;
289 ++t_local_data;
290 }
291 }
292 }
293
294 const Problem *prb_ptr;
295 CHKERR m_field.get_problem(problem, &prb_ptr);
296 boost::shared_ptr<NumeredEntFiniteElement_multiIndex> numered_fes(
298
299 Range in_tets;
300 in_tets.insert_list(data_ptr->evalPointEntityHandle.begin(),
301 data_ptr->evalPointEntityHandle.end());
302 in_tets = in_tets.subset_by_dimension(D);
303
304 if (verb >= NOISY)
305 MOFEM_LOG("SELF", Sev::noisy) << "in tets: " << in_tets << endl;
306
307 auto fe_miit =
308 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().find(
309 finite_element);
310 if (fe_miit !=
311 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().end()) {
312
313 for (auto peit = in_tets.pair_begin(); peit != in_tets.pair_end(); ++peit) {
314
315 auto lo =
316 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().lower_bound(
318 peit->first, (*fe_miit)->getFEUId()));
319 auto hi =
320 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().upper_bound(
322 peit->second, (*fe_miit)->getFEUId()));
323 numered_fes->insert(lo, hi);
324
325 if (verb >= VERY_NOISY)
326 std::cerr << "numered elements:" << std::endl;
327 for (; lo != hi; ++lo)
328 if (verb >= VERY_NOISY)
329 std::cerr << **lo << endl;
330 }
331 if (verb >= VERY_NOISY)
332 std::cerr << std::endl;
333
334 if (auto fe_ptr = data_ptr->feMethodPtr.lock()) {
335
336 MOFEM_LOG_C("FieldEvaluatorSync", Sev::verbose,
337 "Number elements %d to evaluate at proc %d",
338 numered_fes->size(), m_field.get_comm_rank());
340
341 if (!cache_ptr) {
342 cache_ptr = boost::make_shared<CacheTuple>();
343 CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
344
345 MOFEM_LOG("FieldEvaluatorSync", Sev::noisy)
346 << "If you call that function many times in the loop consider to "
347 "set "
348 "cache_ptr outside of the loop. Otherwise code can be slow.";
349 }
350
351 CHKERR m_field.loop_finite_elements(prb_ptr, finite_element, *fe_ptr,
352 lower_rank, upper_rank, numered_fes,
353 bh, cache_ptr, verb);
354
355 } else
356 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
357 "Pointer to element does not exists");
358 }
359
361}
362
364 const double *const point, const double distance, const std::string problem,
365 const std::string finite_element, boost::shared_ptr<SetPtsData> data_ptr,
366 int lower_rank, int upper_rank, boost::shared_ptr<CacheTuple> cache_ptr,
367 MoFEMTypes bh, VERBOSITY_LEVELS verb) {
368 return evalFEAtThePoint<3>(point, distance, problem, finite_element, data_ptr,
369 lower_rank, upper_rank, cache_ptr, bh, verb);
370}
371
373 const double *const point, const double distance, const std::string problem,
374 const std::string finite_element, boost::shared_ptr<SetPtsData> data_ptr,
375 int lower_rank, int upper_rank, boost::shared_ptr<CacheTuple> cache_ptr,
376 MoFEMTypes bh, VERBOSITY_LEVELS verb) {
377 return evalFEAtThePoint<2>(point, distance, problem, finite_element, data_ptr,
378 lower_rank, upper_rank, cache_ptr, bh, verb);
379}
380
381} // namespace MoFEM
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
Definition: LogManager.hpp:345
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
static const double eps
VERBOSITY_LEVELS
Verbosity levels.
Definition: definitions.h:206
@ NOISY
Definition: definitions.h:211
@ VERY_NOISY
Definition: definitions.h:212
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:97
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
FTensor::Index< 'n', SPACE_DIM > n
multi_index_container< boost::shared_ptr< NumeredEntFiniteElement >, indexed_by< ordered_unique< tag< Unique_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_EntFiniteElement, UId, &NumeredEntFiniteElement::getLocalUniqueId > >, ordered_non_unique< tag< Part_mi_tag >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > >, ordered_non_unique< tag< Ent_mi_tag >, const_mem_fun< NumeredEntFiniteElement::interface_type_RefEntity, EntityHandle, &NumeredEntFiniteElement::getEnt > >, ordered_non_unique< tag< Composite_Name_And_Part_mi_tag >, composite_key< NumeredEntFiniteElement, const_mem_fun< NumeredEntFiniteElement::interface_type_FiniteElement, boost::string_ref, &NumeredEntFiniteElement::getNameRef >, member< NumeredEntFiniteElement, unsigned int, &NumeredEntFiniteElement::part > > > > > NumeredEntFiniteElement_multiIndex
MultiIndex for entities for NumeredEntFiniteElement.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
virtual const FiniteElement_multiIndex * get_finite_elements() const =0
Get the finite elements object.
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
FTensor::Index< 'i', SPACE_DIM > i
double D
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
Definition: Templates.hpp:1634
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Core (interface) class.
Definition: Core.hpp:82
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
boost::weak_ptr< SetPtsData > dataPtr
Field evaluator interface.
MoFEMErrorCode buildTree3D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode buildTree2D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode evalFEAtThePoint2D(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > 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 artbitray position.
FieldEvaluatorInterface(const MoFEM::Core &core)
MoFEMErrorCode evalFEAtThePoint3D(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > 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 artbitray position.
MoFEMErrorCode buildTree(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode evalFEAtThePoint(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > 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 artbitray position.
structure to get information form mofem into EntitiesFieldData
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Definition: LogManager.cpp:348
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
Definition: LogManager.cpp:399
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
keeps basic data about problem
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
static MoFEMErrorCode getLocalCoordinatesOnReferenceFourNodeTet(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the Local Coordinates On Reference Four Node Tet object.
Definition: Tools.cpp:90
static MoFEMErrorCode getLocalCoordinatesOnReferenceTriNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the Local Coordinates On Reference Four Node Tet object.
Definition: Tools.cpp:142
base class for all interface classes