12 auto core_log = logging::core::get();
15 "FieldEvaluatorWorld"));
17 "FieldEvaluatorSync"));
19 "FieldEvaluatorSelf"));
30 MOFEM_LOG(
"FieldEvaluatorWorld", Sev::noisy) <<
"Field evaluator intreface";
42 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
49 CHKERR m_field.
get_moab().get_entities_by_dimension(fe_meshset,
D, entities,
55 if (spd_ptr->treePtr) {
56 spd_ptr->treePtr->reset_tree();
58 spd_ptr->treePtr = boost::make_shared<AdaptiveKDTree>(&m_field.
get_moab());
60 CHKERR spd_ptr->treePtr->build_tree(entities, &spd_ptr->rooTreeSet);
66 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
68 return buildTreeImpl<2>(spd_ptr, finite_element,
bit, mask);
73 boost::shared_ptr<SetIntegrationPtsMethodData> spd_ptr,
75 return buildTreeImpl<3>(spd_ptr, finite_element,
bit, mask);
79 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr)
80 : dataPtr(data_ptr) {}
87 if (
auto data_ptr = dataPtr.lock()) {
89 const auto nb_eval_points = data_ptr->nbEvalPoints;
90 const auto &shape_functions = data_ptr->shapeFunctions;
91 const auto &eval_pointentity_handle = data_ptr->evalPointEntityHandle;
93 if (
auto fe_ptr = data_ptr->feMethodPtr) {
98 if (fe_ptr.get() != fe_raw_ptr)
102 const auto fe_ent = fe.numeredEntFiniteElementPtr->getEnt();
104 const auto fe_dim = moab::CN::Dimension(fe_type);
106 auto &gauss_pts = fe.gaussPts;
110 gauss_pts.resize(4, nb_eval_points,
false);
112 &shape_functions(0, 0), &shape_functions(0, 1),
113 &shape_functions(0, 2), &shape_functions(0, 3)};
115 &gauss_pts(0, 0), &gauss_pts(1, 0), &gauss_pts(2, 0)};
118 for (
int nn = 0; nn != nb_eval_points; ++nn) {
119 if (eval_pointentity_handle[nn] == fe_ent) {
120 for (
const int i : {0, 1, 2}) {
121 t_gauss_pts(
i) = t_shape(
i + 1);
123 gauss_pts(3, nb_gauss_pts) = nn;
129 gauss_pts.resize(4, nb_gauss_pts,
true);
130 }
else if (fe_dim == 2) {
131 gauss_pts.resize(3, nb_eval_points,
false);
133 &shape_functions(0, 0), &shape_functions(0, 1),
134 &shape_functions(0, 2)};
136 &gauss_pts(0, 0), &gauss_pts(1, 0)};
138 for (
int nn = 0; nn != nb_eval_points; ++nn) {
139 if (eval_pointentity_handle[nn] == fe_ent) {
140 for (
const int i : {0, 1}) {
141 t_gauss_pts(
i) = t_shape(
i + 1);
143 gauss_pts(2, nb_gauss_pts) = nn;
149 gauss_pts.resize(3, nb_gauss_pts,
true);
152 "Dimension not implemented");
157 <<
"nbEvalOPoints / nbGau_sspt_s: " << nb_eval_points <<
" / "
159 MOFEM_LOG(
"SELF", Sev::noisy) <<
"gauss pts: " << gauss_pts;
164 "Pointer to element does not exists");
168 "Pointer to data does not exists");
175 const double *
const point,
const double distance,
const std::string problem,
176 const std::string finite_element,
177 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
179 boost::shared_ptr<CacheTuple> cache_ptr,
MoFEMTypes bh,
184 std::vector<EntityHandle> leafs_out;
185 CHKERR data_ptr->treePtr->distance_search(point, distance, leafs_out);
187 for (
auto lit : leafs_out)
192 MOFEM_LOG(
"SELF", Sev::noisy) <<
"tree entities: " << tree_ents;
194 data_ptr->evalPointEntityHandle.resize(data_ptr->nbEvalPoints);
195 std::fill(data_ptr->evalPointEntityHandle.begin(),
196 data_ptr->evalPointEntityHandle.end(), 0);
198 std::vector<double> local_coords;
199 std::vector<double> shape;
201 auto nb_eval_points = data_ptr->nbEvalPoints;
203 for (
auto tet : tree_ents) {
207 "Wrong element type, Field Evaluator not implemented for Quads "
212 CHKERR m_field.
get_moab().get_connectivity(tet, conn, num_nodes,
true);
214 if constexpr (
D == 3) {
216 local_coords.resize(3 * nb_eval_points);
217 shape.resize(4 * nb_eval_points);
219 std::array<double, 12> coords;
220 CHKERR m_field.
get_moab().get_coords(conn, num_nodes, coords.data());
223 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
225 CHKERR Tools::shapeFunMBTET<3>(&shape[0], &local_coords[0],
226 &local_coords[1], &local_coords[2],
231 &shape[0], &shape[1], &shape[2], &shape[3]};
233 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
234 &data_ptr->shapeFunctions(0, 2), &data_ptr->shapeFunctions(0, 3)};
238 &local_coords[0], &local_coords[1], &local_coords[2]};
240 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1)),
241 &(data_ptr->localCoords(0, 2))};
243 for (
int n = 0;
n != nb_eval_points; ++
n) {
245 const double eps = data_ptr->eps;
246 if (t_shape(0) >= 0 -
eps && t_shape(0) <= 1 +
eps &&
248 t_shape(1) >= 0 -
eps && t_shape(1) <= 1 +
eps &&
250 t_shape(2) >= 0 -
eps && t_shape(2) <= 1 +
eps &&
252 t_shape(3) >= 0 -
eps && t_shape(3) <= 1 +
eps) {
254 data_ptr->evalPointEntityHandle[
n] = tet;
255 t_shape_data(i4) = t_shape(i4);
256 t_local_data(j3) = t_local(j3);
266 if constexpr (
D == 2) {
268 local_coords.resize(2 * nb_eval_points);
269 shape.resize(3 * nb_eval_points);
271 std::array<double, 9> coords;
272 CHKERR m_field.
get_moab().get_coords(conn, num_nodes, coords.data());
275 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
277 CHKERR Tools::shapeFunMBTRI<2>(&shape[0], &local_coords[0],
278 &local_coords[1], nb_eval_points);
282 &shape[0], &shape[1], &shape[2]};
284 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
285 &data_ptr->shapeFunctions(0, 2)};
289 &local_coords[0], &local_coords[1]};
290 data_ptr->localCoords.resize(nb_eval_points, 2);
292 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1))};
294 for (
int n = 0;
n != nb_eval_points; ++
n) {
296 const double eps = data_ptr->eps;
297 if (t_shape(0) >= 0 -
eps && t_shape(0) <= 1 +
eps &&
299 t_shape(1) >= 0 -
eps && t_shape(1) <= 1 +
eps &&
301 t_shape(2) >= 0 -
eps && t_shape(2) <= 1 +
eps) {
303 data_ptr->evalPointEntityHandle[
n] = tet;
304 t_shape_data(i3) = t_shape(i3);
305 t_local_data(j2) = t_local(j2);
318 boost::shared_ptr<NumeredEntFiniteElement_multiIndex> numered_fes(
322 in_tets.insert_list(data_ptr->evalPointEntityHandle.begin(),
323 data_ptr->evalPointEntityHandle.end());
324 in_tets = in_tets.subset_by_dimension(
D);
327 MOFEM_LOG(
"SELF", Sev::noisy) <<
"in tets: " << in_tets << endl;
335 for (
auto peit = in_tets.pair_begin(); peit != in_tets.pair_end(); ++peit) {
340 peit->first, (*fe_miit)->getFEUId()));
344 peit->second, (*fe_miit)->getFEUId()));
348 MOFEM_LOG(
"FieldEvaluatorSelf", Sev::noisy);
349 for (
auto it = lo; it != hi; ++it) {
350 auto fe_bit = (*it)->getBitRefLevel();
351 if ((fe_bit & mask) != fe_bit)
353 if ((fe_bit &
bit).none())
355 numered_fes->insert(*it);
357 MOFEM_LOG(
"FieldEvaluatorSelf", Sev::noisy) << **it << endl;
360 numered_fes->insert(lo, hi);
362 MOFEM_LOG(
"FieldEvaluatorSelf", Sev::noisy)
363 <<
"numbered elements:" << std::endl;
364 for (; lo != hi; ++lo)
365 MOFEM_LOG(
"FieldEvaluatorSelf", Sev::noisy) << **lo << endl;
370 if (
auto fe_ptr = data_ptr->feMethodPtr) {
373 "Number elements %d to evaluate at proc %d",
377 cache_ptr = boost::make_shared<CacheTuple>();
380 MOFEM_LOG(
"FieldEvaluatorSelf", Sev::noisy)
381 <<
"If you call that function many times in the loop consider to "
382 "set cache_ptr outside of the loop. Otherwise code can be slow.";
386 lower_rank, upper_rank, numered_fes,
387 bh, cache_ptr, verb);
391 "Pointer to element does not exists");
399 const double *
const point,
const double distance,
const std::string problem,
400 const std::string finite_element,
401 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
403 boost::shared_ptr<CacheTuple> cache_ptr,
MoFEMTypes bh,
405 return evalFEAtThePointImpl<2>(point, distance, problem, finite_element,
406 data_ptr, lower_rank, upper_rank,
bit, mask,
407 cache_ptr, bh, verb);
412 const double *
const point,
const double distance,
const std::string problem,
413 const std::string finite_element,
414 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
416 boost::shared_ptr<CacheTuple> cache_ptr,
MoFEMTypes bh,
418 return evalFEAtThePointImpl<3>(point, distance, problem, finite_element,
419 data_ptr, lower_rank, upper_rank,
bit, mask,
420 cache_ptr, bh, verb);
427 const std::string finite_element,
428 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
435 struct Op :
public BaseOP {
439 const std::string finite_element,
440 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
445 : BaseOP(
NOSPACE, BaseOP::OPSPACE),
446 matDataPts(mat_data_pts), finiteElement(finite_element),
447 dataPtr(data_ptr), lowerRank(lower_rank), upperRank(upper_rank),
448 bitRefLevel(
bit), maskRefLevel(mask),
eps(
eps), bhType(bh),
461 auto nb_integration_pts = getGaussPts().size2();
463 auto prb_name = getPtrFE()->problemPtr->getName();
464 auto cache_ptr = getPtrFE()->getCacheWeakPtr();
466 auto coords_at_pts = getCoordsAtGaussPts();
467 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
468 double *coords_ptr = &coords_at_pts(gg, 0);
470 dataPtr->setEvalPoints(coords_ptr, 1);
471 CHKERR field_eval_ptr->evalFEAtThePoint<
D>(
472 coords_ptr,
eps, prb_name, finiteElement, dataPtr, lowerRank,
473 upperRank, bitRefLevel, maskRefLevel, cache_ptr.lock(),
MF_EXIST,
478 auto loop_size = dataPtr->feMethodPtr->getLoopSize();
479 if (loop_size != 1) {
480 if (verbosity >=
QUIET) {
482 LogManager::BitLineID | LogManager::BitScope);
484 <<
"Field evaluator wrong number of elements found: "
489 for (
auto &
m : matDataPts) {
490 auto nb_rows =
m.first->size1();
491 m.second->resize(nb_rows, nb_integration_pts,
false);
494 for (
auto ll = 0; ll != loop_size; ++ll) {
495 for (
int dd = 0; dd != nb_rows; ++dd) {
496 (*
m.second)(dd, gg) += (*
m.first)(dd, ll) / loop_size;
509 std::string finiteElement;
510 boost::shared_ptr<SetIntegrationPtsMethodData> dataPtr;
523 CHKERR data_ptr->treePtr->get_info(data_ptr->rooTreeSet, min, max, dep);
524 double dist = std::sqrt(pow(max[0] - min[0], 2) + pow(max[1] - min[1], 2) +
525 pow(max[2] - min[2], 2));
527 return new Op( mat_data_pts, finite_element, data_ptr,
528 lower_rank, upper_rank,
bit, mask, 1e-6 * dist, bh, verb);
533 boost::shared_ptr<MoFEM::ForcesAndSourcesCore> fe_method_ptr,
534 const double *eval_points,
const int nb_eval_points,
const double eps,
536 : feMethodPtr(fe_method_ptr), evalPoints(eval_points),
537 nbEvalPoints(nb_eval_points),
eps(
eps), verb(verb) {
543 const double *ptr,
const int nb_eval_points) {
545 nbEvalPoints = nb_eval_points;
546 localCoords.resize(nbEvalPoints, 3,
false);
547 shapeFunctions.resize(nbEvalPoints, 4,
false);
552FieldEvaluatorInterface::getDataOperator<2>(
554 const std::string finite_element,
555 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
558 return getDataOperatorImpl<2>(mat_data_pts, finite_element,
559 data_ptr, lower_rank, upper_rank,
bit, mask, bh,
565FieldEvaluatorInterface::getDataOperator<3>(
567 const std::string finite_element,
568 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr,
int lower_rank,
571 return getDataOperatorImpl<3>(mat_data_pts, finite_element,
572 data_ptr, lower_rank, upper_rank,
bit, mask, bh,
#define MOFEM_LOG_C(channel, severity, format,...)
VERBOSITY_LEVELS
Verbosity levels.
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
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
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
auto type_from_handle(const EntityHandle h)
get type from entity handle
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
virtual int get_comm_rank() const =0
Deprecated interface functions.
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble shapeFunctions
void setEvalPoints(const double *ptr, const int nb_eval_points)
SetIntegrationPtsMethodData()=delete
SetIntegrationPtsMethod()=delete
MoFEMErrorCode operator()(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)
Field evaluator interface.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
std::vector< std::pair< boost::shared_ptr< MatrixDouble >, boost::shared_ptr< MatrixDouble > > > MatDataPts
MoFEMErrorCode buildTreeImpl(boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
FieldEvaluatorInterface(const MoFEM::Core &core)
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)
ForcesAndSourcesCore::UserDataOperator * getDataOperatorImpl(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)
structure to get information from mofem into EntitiesFieldData
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
keeps basic data about problem
boost::shared_ptr< NumeredEntFiniteElement_multiIndex > numeredFiniteElementsPtr
store finite elements
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.