18 : cOre(const_cast<
MoFEM::
Core &>(core)), dEbug(false) {}
30 "field < %s > not found, (top tip: check spelling)",
34 auto wrap_lambda_on_enties =
35 [
lambda](boost::shared_ptr<FieldEntity> field_entity_ptr) {
38 auto field_data = field_entity_ptr->getEntFieldData();
39 for (
auto &
v : field_data)
62 "field < %s > not found, (top tip: check spelling)",
67 const auto bit_number = (*fit)->getBitNumber();
76 auto eit = field_ents->get<
Unique_mi_tag>().lower_bound(lo_uid);
78 auto eit_hi = field_ents->get<
Unique_mi_tag>().upper_bound(hi_uid);
80 for (; eit != eit_hi; ++eit) {
88 for (
auto p = ents_ptr->const_pair_begin();
p != ents_ptr->const_pair_end();
90 CHKERR loop_ents(
p->first,
p->second);
95 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
96 get_id_for_max_type<MBENTITYSET>());
104 const std::string field_name_y,
bool error_if_missing,
Range *ents_ptr) {
112 "x field < %s > not found, (top tip: check spelling)",
113 field_name_x.c_str());
118 "y field < %s > not found, (top tip: check spelling)",
119 field_name_y.c_str());
121 if ((*x_fit)->getSpace() != (*y_fit)->getSpace()) {
123 "space for field < %s > and field <%s> are not compatible for "
125 field_name_x.c_str(), field_name_y.c_str());
127 if ((*x_fit)->getNbOfCoeffs() != (*y_fit)->getNbOfCoeffs()) {
130 "rank for field < %s > and field <%s> are not compatible for fieldblas",
131 field_name_x.c_str(), field_name_y.c_str());
134 auto wrap_lambda_on_enties =
135 [
lambda](boost::shared_ptr<FieldEntity> y_field_entity_ptr,
136 const boost::shared_ptr<FieldEntity> x_field_entity_ptr) {
139 auto x_field_data = x_field_entity_ptr->getEntFieldData();
140 auto y_field_data = y_field_entity_ptr->getEntFieldData();
141 const auto size_x = x_field_data.size();
142 const auto size_y = y_field_data.size();
145 for (; dd != std::min(size_x, size_y); ++dd)
146 y_field_data[dd] =
lambda(y_field_data[dd], x_field_data[dd]);
147 for (; dd < size_y; ++dd)
148 y_field_data[dd] = 0;
154 field_name_y, error_if_missing, ents_ptr);
160 const std::string field_name_x,
161 const std::string field_name_y,
162 bool error_if_missing,
Range *ents_ptr) {
172 "x field < %s > not found, (top tip: check spelling)",
173 field_name_x.c_str());
178 "y field < %s > not found, (top tip: check spelling)",
179 field_name_y.c_str());
183 const auto x_bit_number = (*x_fit)->getBitNumber();
184 const auto y_bit_number = (*y_fit)->getBitNumber();
185 const auto y_eit_hi = field_ents->get<
Unique_mi_tag>().upper_bound(
191 const auto x_lo_uid =
193 const auto x_hi_uid =
197 auto x_eit = field_ents->get<
Unique_mi_tag>().lower_bound(x_lo_uid);
199 auto x_eit_hi = field_ents->get<
Unique_mi_tag>().upper_bound(x_hi_uid);
201 for (; x_eit != x_eit_hi;) {
204 (*y_fit)->getBitNumber(), (*x_eit)->getEnt());
205 auto y_eit = field_ents->get<
Unique_mi_tag>().find(y_lo_uid);
207 if (y_eit == field_ents->end()) {
209 if (error_if_missing) {
211 "Missing entity in y_field.");
219 if (x_eit == x_eit_hi)
221 if (y_eit == y_eit_hi)
223 if ((*y_eit)->getEnt() != (*x_eit)->getEnt())
243 for (
auto p = ents_ptr->const_pair_begin();
p != ents_ptr->const_pair_end();
245 CHKERR loop_ents(
p->first,
p->second);
250 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
251 get_id_for_max_type<MBENTITYSET>());
258 const std::string field_name_x,
259 const std::string field_name_y,
260 bool error_if_missing,
261 bool create_if_missing) {
262 if (create_if_missing)
264 <<
"Option create_if_missing is set to true but have no effect";
270 const std::string field_name_x,
271 const std::string field_name_y,
272 bool error_if_missing,
Range *ents_ptr) {
274 auto axpy = [alpha](
const double fy,
const double fx) {
275 return fy + alpha * fx;
283 const std::string field_name_x,
284 const std::string field_name_y,
285 bool error_if_missing,
286 bool create_if_missing) {
287 return fieldAxpy(alpha, field_name_x, field_name_y, error_if_missing);
291 const std::string field_name_x,
292 const std::string field_name_y,
293 bool error_if_missing,
Range *ents_ptr) {
295 auto copy = [alpha](
const double fy,
const double fx) {
return alpha * fx; };
302 const std::string field_name_x,
303 const std::string field_name_y,
304 bool error_if_missing,
305 bool create_if_missing) {
306 return fieldCopy(alpha, field_name_x, field_name_y, error_if_missing);
314 auto scale = [alpha](
const double v) {
return v * alpha; };
328 CHKERR m_field.
get_moab().get_entities_by_type(meshset, MBVERTEX, verts,
331 verts = intersect(*sub_verts, verts);
337 count(0), total(0) {}
344 if (*vit != entPtr->getEnt())
346 "Inconsistent entity %ld != %ld", *vit, entPtr->getEnt());
348 CHKERR mField.get_moab().coords_iterate(vit, verts.end(), xCoord,
349 yCoord, zCoord, count);
350 CHKERR lambda(entPtr->getEntFieldData(), xCoord, yCoord, zCoord);
361 if (total != verts.size())
363 "Inconsistent number of vertices in field meshset and in the "
365 total, verts.size());
400 auto dit = dofs_ptr->get<
Unique_mi_tag>().lower_bound(lo_uid);
401 auto hi_dit = dofs_ptr->get<
Unique_mi_tag>().upper_bound(hi_uid);
403 for (; dit != hi_dit; dit++)
404 (*dit)->getFieldData() = val;
422 auto dit = dofs_ptr->get<
Unique_mi_tag>().lower_bound(lo_uid);
423 auto hi_dit = dofs_ptr->get<
Unique_mi_tag>().upper_bound(hi_uid);
427 for (; dit != hi_dit; dit++) {
428 ent = (*dit)->getEnt();
430 if (ents.find(ent) == ents.end()) {
438 (*dit)->getFieldData() = val;
452 " field < %s > not found, (top tip: check spelling)",
460 for (; dit != hi_dit; dit++) {
461 (*dit)->getFieldData() = val;
#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 ...
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
#define MOFEM_LOG(channel, severity)
Log.
const double v
phase velocity of light in medium (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
EntityHandle get_id_for_max_type()
EntityHandle get_id_for_min_type()
constexpr auto field_name
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
virtual moab::Interface & get_moab()=0
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
Deprecated interface functions.
Data structure to exchange data between mofem and User Loop Methods on entities.
DEPRECATED MoFEMErrorCode fieldLambda(TwoFieldFunctionOnValues lambda, const std::string field_name_x, const std::string field_name_y, bool error_if_missing, bool create_if_missing)
MoFEMErrorCode fieldAxpy(const double alpha, const std::string field_name_x, const std::string field_name_y, bool error_if_missing=false, Range *ents_ptr=nullptr)
axpy fields
boost::function< MoFEMErrorCode(boost::shared_ptr< FieldEntity > field_entity_ptr)> OneFieldFunctionOnEntities
MoFEMErrorCode fieldLambdaOnEntities(OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
boost::function< double(const double y_field_value_reference, const double x_field_value)> TwoFieldFunctionOnValues
MoFEMErrorCode setVertexDofs(VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr)
Set DOFs on vertices using user function.
boost::function< MoFEMErrorCode(VectorAdaptor &&field_data, double *xcoord, double *ycoord, double *zcoord)> VertexCoordsFunction
MoFEMErrorCode setField(const double val, const EntityType type, const std::string field_name)
boost::function< double(const double field_val)> OneFieldFunctionOnValues
function to set a field value
MoFEMErrorCode fieldScale(const double alpha, const std::string field_name, Range *ents_ptr=nullptr)
field scale
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
FieldBlas(const MoFEM::Core &core)
const MoFEM::Interface & cOre
boost::function< MoFEMErrorCode(boost::shared_ptr< FieldEntity > y_field_entity_ptr, const boost::shared_ptr< FieldEntity > x_field_entity_ptr)> TwoFieldFunctionOnEntities
MoFEMErrorCode fieldCopy(const double alpha, const std::string field_name_x, const std::string field_name_y, bool error_if_missing=false, Range *ents_ptr=nullptr)
copy and scale fields
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
static UId getLoLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
static UId getHiLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
MultiIndex Tag for field name.
base class for all interface classes