v0.15.0
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
MoFEM::FieldEvaluatorInterface Struct Reference

Field evaluator interface. More...

#include "src/interfaces/FieldEvaluator.hpp"

Inheritance diagram for MoFEM::FieldEvaluatorInterface:
[legend]
Collaboration diagram for MoFEM::FieldEvaluatorInterface:
[legend]

Classes

struct  SetIntegrationPtsMethod
 
struct  SetIntegrationPtsMethodData
 

Public Types

using SetPtsData = SetIntegrationPtsMethodData
 
using MatDataPts = std::vector< std::pair< boost::shared_ptr< MatrixDouble >, boost::shared_ptr< MatrixDouble > > >
 

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FieldEvaluatorInterface (const MoFEM::Core &core)
 
template<typename VE , typename SPD = SetIntegrationPtsMethodData, typename SP = SetIntegrationPtsMethod>
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.
 
template<int D>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element)
 Build spatial tree.
 
template<int D>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 Build spatial tree.
 
template<int D>
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 arbitrary position.
 
template<int D>
ForcesAndSourcesCore::UserDataOperatorgetDataOperator (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.
 
template<int D>
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 arbitrary position.
 
template<>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 
template<>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 
template<>
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, VERBOSITY_LEVELS verb)
 
template<>
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, VERBOSITY_LEVELS verb)
 
template<>
ForcesAndSourcesCore::UserDataOperatorgetDataOperator (FieldEvaluatorInterface::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, VERBOSITY_LEVELS verb)
 
template<>
ForcesAndSourcesCore::UserDataOperatorgetDataOperator (FieldEvaluatorInterface::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, VERBOSITY_LEVELS verb)
 
template<>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 
template<>
MoFEMErrorCode buildTree (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 
template<>
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, VERBOSITY_LEVELS verb)
 
template<>
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, VERBOSITY_LEVELS verb)
 
template<>
ForcesAndSourcesCore::UserDataOperatorgetDataOperator (FieldEvaluatorInterface::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, VERBOSITY_LEVELS verb)
 
template<>
ForcesAndSourcesCore::UserDataOperatorgetDataOperator (FieldEvaluatorInterface::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, VERBOSITY_LEVELS verb)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 

Private Member Functions

template<int D>
MoFEMErrorCode buildTreeImpl (boost::shared_ptr< SetIntegrationPtsMethodData > spd_ptr, const std::string finite_element, BitRefLevel bit, BitRefLevel mask)
 
template<int D>
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)
 
template<int D>
ForcesAndSourcesCore::UserDataOperatorgetDataOperatorImpl (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)
 

Private Attributes

MoFEM::CorecOre
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 

Detailed Description

Field evaluator interface.

Examples
mofem/atom_tests/field_evaluator.cpp, mofem/tutorials/adv-0/plastic.cpp, mofem/tutorials/adv-2/thermo_elastic.cpp, mofem/tutorials/adv-4/dynamic_first_order_con_law.cpp, mofem/tutorials/adv-5/seepage.cpp, mofem/tutorials/adv-6/between_meshes_dg_projection.cpp, mofem/tutorials/max-1/lorentz_force.cpp, mofem/tutorials/vec-7/adjoint.cpp, mofem/users_modules/basic_finite_elements/elasticity/elasticity.cpp, plastic.cpp, thermo_elastic.cpp, and thermoplastic.cpp.

Definition at line 21 of file FieldEvaluator.hpp.

Member Typedef Documentation

◆ MatDataPts

using MoFEM::FieldEvaluatorInterface::MatDataPts = std::vector< std::pair<boost::shared_ptr<MatrixDouble>, boost::shared_ptr<MatrixDouble> > >

Definition at line 145 of file FieldEvaluator.hpp.

◆ SetPtsData

Definition at line 24 of file FieldEvaluator.hpp.

Constructor & Destructor Documentation

◆ FieldEvaluatorInterface()

MoFEM::FieldEvaluatorInterface::FieldEvaluatorInterface ( const MoFEM::Core core)

Definition at line 8 of file FieldEvaluator.cpp.

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}
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.
Core (interface) class.
Definition Core.hpp:82
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.

Member Function Documentation

◆ buildTree() [1/6]

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element 
)
inline

Build spatial tree.

Parameters
finite_elementfinite element name
Returns
MoFEMErrorCode
Examples
mofem/atom_tests/field_evaluator.cpp, and mofem/tutorials/max-1/lorentz_force.cpp.

Definition at line 70 of file FieldEvaluator.hpp.

71 {
72 return buildTree<D>(spd_ptr, finite_element, BitRefLevel(), BitRefLevel());
73 }
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40

◆ buildTree() [2/6]

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)

Build spatial tree.

Parameters
finite_elementfinite element name
bitbit reference level
maskmask reference level
Returns
MoFEMErrorCode

◆ buildTree() [3/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)

◆ buildTree() [4/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)

◆ buildTree() [5/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)

Definition at line 61 of file FieldEvaluator.cpp.

63 {
64 return buildTreeImpl<2>(spd_ptr, finite_element, bit, mask);
65}
auto bit
set bit

◆ buildTree() [6/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTree ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)

Definition at line 68 of file FieldEvaluator.cpp.

70 {
71 return buildTreeImpl<3>(spd_ptr, finite_element, bit, mask);
72}

◆ buildTreeImpl()

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::buildTreeImpl ( boost::shared_ptr< SetIntegrationPtsMethodData spd_ptr,
const std::string  finite_element,
BitRefLevel  bit,
BitRefLevel  mask 
)
private

Definition at line 41 of file FieldEvaluator.cpp.

43 {
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 if (bit.any()) {
52 CHKERR m_field.getInterface<BitRefManager>()->filterEntitiesByRefLevel(
53 bit, mask, entities);
54 }
55 spd_ptr->treePtr.reset(new AdaptiveKDTree(&m_field.get_moab()));
56 CHKERR spd_ptr->treePtr->build_tree(entities, &spd_ptr->rooTreeSet);
58}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
double D

◆ evalFEAtThePoint() [1/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

◆ evalFEAtThePoint() [2/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

◆ evalFEAtThePoint() [3/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

Definition at line 394 of file FieldEvaluator.cpp.

400 {
401 return evalFEAtThePointImpl<2>(point, distance, problem, finite_element,
402 data_ptr, lower_rank, upper_rank, bit, mask,
403 cache_ptr, bh, verb);
404}

◆ evalFEAtThePoint() [4/6]

template<>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

Definition at line 407 of file FieldEvaluator.cpp.

413 {
414 return evalFEAtThePointImpl<3>(point, distance, problem, finite_element,
415 data_ptr, lower_rank, upper_rank, bit, mask,
416 cache_ptr, bh, verb);
417}

◆ evalFEAtThePoint() [5/6]

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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 arbitrary position.

Parameters
pointpoint used to find tetrahedrons
distancedistance from the point where tetrahedrons are searched
problemproblem name
finite_elementfinite element name
data_ptrpointer to data abut gauss points
lower_ranklower processor
upper_rankupper process
bitbit reference level
maskmask reference level
cache_ptrcache or problem entities
bhcontrol looping over entities, e.g. throwing error if element not found
verbverbosity level
Returns
MoFEMErrorCode

◆ evalFEAtThePoint() [6/6]

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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 
)
inline

Evaluate field at arbitrary position.

std::array<double, 3> point = {0, 0, 0};
const double dist = 0.3;
std::array<double, 6> eval_points = {-1., -1., -1., 1., 1., 1. };
using VolEle = VolumeElementForcesAndSourcesCore;
auto data_ptr =
m_field.getInterface<FieldEvaluatorInterface>()->
getData<VolEle>(point.data(), point.size/3);
if(auto vol_ele = data_ptr->feMethod.lock()) {
// push operators to finite element instance, e.g.
vol_ele->getOpPtrVector().push_back(new MyOp());
// iterate over elemnts with evaluated points
auto cache_ptr = boost::make_shared<CacheTuple>();
CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
// SPACE_DIM is 3 in this example
CHKERR m_field.getInterface<FieldEvaluatorInterface>()
->evalFEAtThePoint<3>(point.data(), dist, prb_ptr->getName(),
"FINITE_ELEMENT_NAME",
data_ptr, m_field.get_comm_rank(),
m_field.get_comm_rank(), cache_ptr);
}
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
Parameters
pointpoint used to find tetrahedrons
distancedistance from the point where tetrahedrons are searched
problemproblem name
finite_elementfinite element name
data_ptrpointer to data abut gauss points
lower_ranklower processor
upper_rankupper process
cache_ptrcache or problem entities
bhcontrol looping over entities, e.g. throwing error if element not found
verbverbosity level
Returns
MoFEMErrorCode
Examples
mofem/atom_tests/field_evaluator.cpp, and mofem/tutorials/max-1/lorentz_force.cpp.

Definition at line 134 of file FieldEvaluator.hpp.

139 {
140 return evalFEAtThePoint<D>(point, distance, problem, finite_element,
141 data_ptr, lower_rank, upper_rank, BitRefLevel(),
142 BitRefLevel(), cache_ptr, bh, verb);
143 }

◆ evalFEAtThePointImpl()

template<int D>
MoFEMErrorCode MoFEM::FieldEvaluatorInterface::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 
)
private

Definition at line 170 of file FieldEvaluator.cpp.

176 {
177 CoreInterface &m_field = cOre;
179
180 std::vector<EntityHandle> leafs_out;
181 CHKERR data_ptr->treePtr->distance_search(point, distance, leafs_out);
182 Range tree_ents;
183 for (auto lit : leafs_out) //{
184 CHKERR m_field.get_moab().get_entities_by_dimension(lit, D, tree_ents,
185 true);
186
187 if (verb >= NOISY)
188 MOFEM_LOG("SELF", Sev::noisy) << "tree entities: " << tree_ents;
189
190 data_ptr->evalPointEntityHandle.resize(data_ptr->nbEvalPoints);
191 std::fill(data_ptr->evalPointEntityHandle.begin(),
192 data_ptr->evalPointEntityHandle.end(), 0);
193
194 std::vector<double> local_coords;
195 std::vector<double> shape;
196
197 auto nb_eval_points = data_ptr->nbEvalPoints;
198
199 for (auto tet : tree_ents) {
200 if ((type_from_handle(tet) != moab::MBTRI) &&
201 (type_from_handle(tet) != moab::MBTET)) {
202 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
203 "Wrong element type, Field Evaluator not implemented for Quads "
204 "and Hexes");
205 }
206 const EntityHandle *conn;
207 int num_nodes;
208 CHKERR m_field.get_moab().get_connectivity(tet, conn, num_nodes, true);
209
210 if constexpr (D == 3) {
211
212 local_coords.resize(3 * nb_eval_points);
213 shape.resize(4 * nb_eval_points);
214
215 std::array<double, 12> coords;
216 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
217
219 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
220 &local_coords[0]);
221 CHKERR Tools::shapeFunMBTET<3>(&shape[0], &local_coords[0],
222 &local_coords[1], &local_coords[2],
223 nb_eval_points);
224
225 FTensor::Index<'i', 4> i4;
227 &shape[0], &shape[1], &shape[2], &shape[3]};
229 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
230 &data_ptr->shapeFunctions(0, 2), &data_ptr->shapeFunctions(0, 3)};
231
232 FTensor::Index<'j', 3> j3;
234 &local_coords[0], &local_coords[1], &local_coords[2]};
236 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1)),
237 &(data_ptr->localCoords(0, 2))};
238
239 for (int n = 0; n != nb_eval_points; ++n) {
240
241 const double eps = data_ptr->eps;
242 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
243
244 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
245
246 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps &&
247
248 t_shape(3) >= 0 - eps && t_shape(3) <= 1 + eps) {
249
250 data_ptr->evalPointEntityHandle[n] = tet;
251 t_shape_data(i4) = t_shape(i4);
252 t_local_data(j3) = t_local(j3);
253 }
254
255 ++t_shape;
256 ++t_shape_data;
257 ++t_local;
258 ++t_local_data;
259 }
260 }
261
262 if constexpr (D == 2) {
263
264 local_coords.resize(2 * nb_eval_points);
265 shape.resize(3 * nb_eval_points);
266
267 std::array<double, 9> coords;
268 CHKERR m_field.get_moab().get_coords(conn, num_nodes, coords.data());
269
271 coords.data(), &data_ptr->evalPoints[0], nb_eval_points,
272 &local_coords[0]);
273 CHKERR Tools::shapeFunMBTRI<2>(&shape[0], &local_coords[0],
274 &local_coords[1], nb_eval_points);
275
276 FTensor::Index<'i', 3> i3;
278 &shape[0], &shape[1], &shape[2]};
280 &data_ptr->shapeFunctions(0, 0), &data_ptr->shapeFunctions(0, 1),
281 &data_ptr->shapeFunctions(0, 2)};
282
283 FTensor::Index<'j', 2> j2;
285 &local_coords[0], &local_coords[1]};
286 data_ptr->localCoords.resize(nb_eval_points, 2);
288 &(data_ptr->localCoords(0, 0)), &(data_ptr->localCoords(0, 1))};
289
290 for (int n = 0; n != nb_eval_points; ++n) {
291
292 const double eps = data_ptr->eps;
293 if (t_shape(0) >= 0 - eps && t_shape(0) <= 1 + eps &&
294
295 t_shape(1) >= 0 - eps && t_shape(1) <= 1 + eps &&
296
297 t_shape(2) >= 0 - eps && t_shape(2) <= 1 + eps) {
298
299 data_ptr->evalPointEntityHandle[n] = tet;
300 t_shape_data(i3) = t_shape(i3);
301 t_local_data(j2) = t_local(j2);
302 }
303
304 ++t_shape;
305 ++t_shape_data;
306 ++t_local;
307 ++t_local_data;
308 }
309 }
310 }
311
312 const Problem *prb_ptr;
313 CHKERR m_field.get_problem(problem, &prb_ptr);
314 boost::shared_ptr<NumeredEntFiniteElement_multiIndex> numered_fes(
316
317 Range in_tets;
318 in_tets.insert_list(data_ptr->evalPointEntityHandle.begin(),
319 data_ptr->evalPointEntityHandle.end());
320 in_tets = in_tets.subset_by_dimension(D);
321
322 if (verb >= NOISY)
323 MOFEM_LOG("SELF", Sev::noisy) << "in tets: " << in_tets << endl;
324
325 auto fe_miit =
326 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().find(
327 finite_element);
328 if (fe_miit !=
329 m_field.get_finite_elements()->get<FiniteElement_name_mi_tag>().end()) {
330
331 for (auto peit = in_tets.pair_begin(); peit != in_tets.pair_end(); ++peit) {
332
333 auto lo =
334 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().lower_bound(
336 peit->first, (*fe_miit)->getFEUId()));
337 auto hi =
338 prb_ptr->numeredFiniteElementsPtr->get<Unique_mi_tag>().upper_bound(
340 peit->second, (*fe_miit)->getFEUId()));
341
342 if (bit.any()) {
343 if (verb >= VERY_NOISY)
344 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy);
345 for (auto it = lo; it != hi; ++it) {
346 auto fe_bit = (*it)->getBitRefLevel();
347 if ((fe_bit & mask) != fe_bit)
348 continue;
349 if ((fe_bit & bit).none())
350 continue;
351 numered_fes->insert(*it);
352 if (verb >= VERY_NOISY)
353 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy) << **it << endl;
354 }
355 } else {
356 numered_fes->insert(lo, hi);
357 if (verb >= VERY_NOISY) {
358 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy)
359 << "numbered elements:" << std::endl;
360 for (; lo != hi; ++lo)
361 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy) << **lo << endl;
362 }
363 }
364 }
365
366 if (auto fe_ptr = data_ptr->feMethodPtr.lock()) {
367
368 MOFEM_LOG_C("FieldEvaluatorSelf", Sev::noisy,
369 "Number elements %d to evaluate at proc %d",
370 numered_fes->size(), m_field.get_comm_rank());
371
372 if (!cache_ptr) {
373 cache_ptr = boost::make_shared<CacheTuple>();
374 CHKERR m_field.cache_problem_entities(prb_ptr->getName(), cache_ptr);
375
376 MOFEM_LOG("FieldEvaluatorSelf", Sev::noisy)
377 << "If you call that function many times in the loop consider to "
378 "set cache_ptr outside of the loop. Otherwise code can be slow.";
379 }
380
381 CHKERR m_field.loop_finite_elements(prb_ptr, finite_element, *fe_ptr,
382 lower_rank, upper_rank, numered_fes,
383 bh, cache_ptr, verb);
384
385 } else
386 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
387 "Pointer to element does not exists");
388 }
389
391}
#define MOFEM_LOG_C(channel, severity, format,...)
static const double eps
@ NOISY
@ VERY_NOISY
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
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.
const double n
refractive index of diffusive medium
auto type_from_handle(const EntityHandle h)
get type from entity handle
UId getLocalUniqueIdCalculate() const
Generate UId for finite element entity.
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:88
static MoFEMErrorCode getLocalCoordinatesOnReferenceThreeNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the local coordinates on reference three node tri object.
Definition Tools.cpp:188

◆ getData()

template<typename VE , typename SPD , typename SP >
boost::shared_ptr< SPD > MoFEM::FieldEvaluatorInterface::getData ( const double ptr = nullptr,
const int  nb_eval_points = 0,
const double  eps = 1e-12,
VERBOSITY_LEVELS  verb = QUIET 
)

Get the Data object.

Pack pointers with data structures for field evaluator and finite element. Function return shared pointer if returned shared pointer is reset; all data are destroyed. The idea is to pack all data in one structure, create shared pointer to it and return aliased shared pointer to one of the elements of the data structure. It is a bit complicated, but has great flexibility.

Template Parameters
VE
SetIntegrationPtsMethodData
SetIntegrationPtsMethod
Parameters
ptr
nb_eval_points
eps
verb
Returns
boost::shared_ptr<SPD>
Examples
mofem/atom_tests/field_evaluator.cpp, mofem/tutorials/adv-6/between_meshes_dg_projection.cpp, and mofem/tutorials/max-1/lorentz_force.cpp.

Definition at line 324 of file FieldEvaluator.hpp.

325 {
326
327 struct PackData {
328 boost::scoped_ptr<VE> elePtr; //< finite element pointer
329 boost::scoped_ptr<SPD> setPtsDataPtr; //< integration pts data
330 boost::scoped_ptr<SP> setPtsPtr; //< integration pts method
331 };
332 boost::shared_ptr<PackData> pack_data(new PackData());
333 MoFEM::Interface &m_field = cOre;
334 pack_data->elePtr.reset(new VE(m_field));
335
336 pack_data->setPtsDataPtr.reset(
337 new SPD(boost::shared_ptr<VE>(pack_data, pack_data->elePtr.get()), ptr,
338 nb_eval_points, eps, verb));
339 pack_data->setPtsPtr.reset(new SP(
340 boost::shared_ptr<SPD>(pack_data, pack_data->setPtsDataPtr.get())));
341 pack_data->elePtr->setRuleHook = boost::ref(*pack_data->setPtsPtr);
342 boost::shared_ptr<SPD> data(pack_data, pack_data->setPtsDataPtr.get());
343 return data;
344}
Deprecated interface functions.

◆ getDataOperator() [1/5]

template<>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperator ( FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

◆ getDataOperator() [2/5]

template<>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperator ( FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

◆ getDataOperator() [3/5]

template<>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperator ( FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

Definition at line 547 of file FieldEvaluator.cpp.

553 {
554 return getDataOperatorImpl<2>(mat_data_pts, finite_element,
555 data_ptr, lower_rank, upper_rank, bit, mask, bh,
556 verb);
557}

◆ getDataOperator() [4/5]

template<>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperator ( FieldEvaluatorInterface::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,
VERBOSITY_LEVELS  verb 
)

Definition at line 560 of file FieldEvaluator.cpp.

566 {
567 return getDataOperatorImpl<3>(mat_data_pts, finite_element,
568 data_ptr, lower_rank, upper_rank, bit, mask, bh,
569 verb);
570}

◆ getDataOperator() [5/5]

template<int D>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperator ( 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.

Parameters
vec_data_ptsvector of data shared between physical and post proc elements
mat_data_ptsmatrix of data shared between physical and post proc elements
finite_elementfinite element name on which field is evaluated
data_ptrdata about field evaluator
bhcontrol looping over entities, e.g. throwing error if element not found

◆ getDataOperatorImpl()

template<int D>
ForcesAndSourcesCore::UserDataOperator * MoFEM::FieldEvaluatorInterface::getDataOperatorImpl ( FieldEvaluatorInterface::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 
)
private

Definition at line 421 of file FieldEvaluator.cpp.

427 {
428
429 using BaseOP = ForcesAndSourcesCore::UserDataOperator;
430
431 struct Op : public BaseOP {
432 Op(
433
434 FieldEvaluatorInterface::MatDataPts mat_data_pts,
435 const std::string finite_element,
436 boost::shared_ptr<SetIntegrationPtsMethodData> data_ptr, int lower_rank,
437 int upper_rank, BitRefLevel bit, BitRefLevel mask, double eps,
439
440 )
441 : BaseOP(NOSPACE, BaseOP::OPSPACE),
442 matDataPts(mat_data_pts), finiteElement(finite_element),
443 dataPtr(data_ptr), lowerRank(lower_rank), upperRank(upper_rank),
444 bitRefLevel(bit), maskRefLevel(mask), eps(eps), bhType(bh),
445 verbosity(verb) {}
446
447 MoFEMErrorCode doWork(int side, EntityType type,
448 EntitiesFieldData::EntData &data) {
450
451 MOFEM_LOG_CHANNEL("SELF");
452 MOFEM_LOG_TAG("SELF", "FieldEvaluator");
453
454 auto &m_field = this->getPtrFE()->mField;
455 auto field_eval_ptr = m_field.getInterface<FieldEvaluatorInterface>();
456
457 auto nb_integration_pts = getGaussPts().size2();
458
459 auto prb_name = getPtrFE()->problemPtr->getName();
460 auto cache_ptr = getPtrFE()->getCacheWeakPtr();
461
462 auto coords_at_pts = getCoordsAtGaussPts();
463 for (int gg = 0; gg != nb_integration_pts; ++gg) {
464 double *coords_ptr = &coords_at_pts(gg, 0);
465
466 dataPtr->setEvalPoints(coords_ptr, 1);
467 CHKERR field_eval_ptr->evalFEAtThePoint<D>(
468 coords_ptr, eps, prb_name, finiteElement, dataPtr, lowerRank,
469 upperRank, bitRefLevel, maskRefLevel, cache_ptr.lock(), MF_EXIST,
470 QUIET);
471
472 // Get number of found finite elements, that consist gauss point.
473 // E.g. for edge flip expected is that only one element will be found.
474 auto loop_size = dataPtr->feMethodPtr.lock()->getLoopSize();
475 if (loop_size != 1) {
476 if (verbosity >= QUIET) {
478 LogManager::BitLineID | LogManager::BitScope);
479 MOFEM_LOG("SELF", Sev::warning)
480 << "Field evaluator wrong number of elements found: "
481 << loop_size;
482 }
483 }
484
485 for (auto &m : matDataPts) {
486 auto nb_rows = m.first->size1();
487 m.second->resize(nb_rows, nb_integration_pts, false);
488 if (gg == 0)
489 m.second->clear();
490 for (auto ll = 0; ll != loop_size; ++ll) {
491 for (int dd = 0; dd != nb_rows; ++dd) {
492 (*m.second)(dd, gg) += (*m.first)(dd, ll) / loop_size;
493 }
494 }
495 }
496 }
497
498 MOFEM_LOG_CHANNEL("SELF");
499
501 }
502
503 private:
504 FieldEvaluatorInterface::MatDataPts matDataPts;
505 std::string finiteElement;
506 boost::shared_ptr<SetIntegrationPtsMethodData> dataPtr;
507 int lowerRank;
508 int upperRank;
509 BitRefLevel bitRefLevel;
510 BitRefLevel maskRefLevel;
511 double eps;
512 MoFEMTypes bhType;
513 VERBOSITY_LEVELS verbosity;
514 };
515
516 double min[3];
517 double max[3];
518 unsigned int dep;
519 CHKERR data_ptr->treePtr->get_info(data_ptr->rooTreeSet, min, max, dep);
520 double dist = std::sqrt(pow(max[0] - min[0], 2) + pow(max[1] - min[1], 2) +
521 pow(max[2] - min[2], 2));
522
523 return new Op( mat_data_pts, finite_element, data_ptr,
524 lower_rank, upper_rank, bit, mask, 1e-6 * dist, bh, verb);
525}
VERBOSITY_LEVELS
Verbosity levels.
@ QUIET
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ MF_EXIST
@ NOSPACE
Definition definitions.h:83
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition ddTensor0.hpp:33
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
FTensor::Index< 'm', 3 > m

◆ query_interface()

MoFEMErrorCode MoFEM::FieldEvaluatorInterface::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Implements MoFEM::UnknownInterface.

Definition at line 33 of file FieldEvaluator.cpp.

34 {
36 *iface = const_cast<FieldEvaluatorInterface *>(this);
38}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FieldEvaluatorInterface(const MoFEM::Core &core)

Member Data Documentation

◆ cOre

MoFEM::Core& MoFEM::FieldEvaluatorInterface::cOre
private

Definition at line 201 of file FieldEvaluator.hpp.


The documentation for this struct was generated from the following files: