|
| v0.14.0
|
Basic algebra on fields.
More...
#include <src/interfaces/FieldBlas.hpp>
|
MoFEMErrorCode | query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const |
|
| FieldBlas (const MoFEM::Core &core) |
|
| ~FieldBlas ()=default |
|
MoFEMErrorCode | fieldLambdaOnEntities (OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr) |
| filed lambda More...
|
|
MoFEMErrorCode | fieldLambdaOnValues (OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr) |
| filed lambda More...
|
|
MoFEMErrorCode | fieldLambdaOnValues (TwoFieldFunctionOnValues lambda, const std::string field_name_x, const std::string field_name_y, bool error_if_missing=false, Range *ents_ptr=nullptr) |
| filed lambda More...
|
|
MoFEMErrorCode | fieldLambdaOnEntities (TwoFieldFunctionOnEntities lambda, const std::string field_name_x, const std::string field_name_y, bool error_if_missing=false, Range *ents_ptr=nullptr) |
| filed lambda More...
|
|
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 More...
|
|
DEPRECATED MoFEMErrorCode | fieldAxpy (const double alpha, const std::string field_name_x, const std::string field_name_y, bool error_if_missing, bool create_if_missing) |
|
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 More...
|
|
DEPRECATED MoFEMErrorCode | fieldCopy (const double alpha, const std::string field_name_x, const std::string field_name_y, bool error_if_missing, bool create_if_missing) |
|
MoFEMErrorCode | fieldScale (const double alpha, const std::string field_name, Range *ents_ptr=nullptr) |
| field scale More...
|
|
MoFEMErrorCode | setVertexDofs (VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr) |
| Set DOFs on vertices using user function. More...
|
|
MoFEMErrorCode | setField (const double val, const EntityType type, const std::string field_name) |
|
MoFEMErrorCode | setField (const double val, const EntityType type, const Range &ents, const std::string field_name) |
| set field More...
|
|
MoFEMErrorCode | setField (const double val, const std::string field_name) |
| set field More...
|
|
template<class IFACE > |
MoFEMErrorCode | registerInterface (bool error_if_registration_failed=true) |
| Register interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE *&iface) const |
| Get interface reference to pointer of interface. More...
|
|
template<class IFACE > |
MoFEMErrorCode | getInterface (IFACE **const iface) const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get interface pointer to pointer of interface. More...
|
|
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0> |
IFACE | getInterface () const |
| Get reference to interface. More...
|
|
template<class IFACE > |
IFACE * | getInterface () const |
| Function returning pointer to interface. More...
|
|
virtual | ~UnknownInterface ()=default |
|
Basic algebra on fields.
- Examples
- adolc_plasticity.cpp, approx_sphere.cpp, dm_partitioned_no_field.cpp, dynamic_first_order_con_law.cpp, elasticity.cpp, field_blas.cpp, free_surface.cpp, heat_equation.cpp, level_set.cpp, mesh_smoothing.cpp, navier_stokes.cpp, nonlinear_dynamics.cpp, reaction_diffusion.cpp, Remodeling.cpp, simple_contact_thermal.cpp, test_jacobian_of_simple_contact_element.cpp, and thermo_elastic.cpp.
Definition at line 21 of file FieldBlas.hpp.
◆ OneFieldFunctionOnEntities
◆ OneFieldFunctionOnValues
◆ TwoFieldFunctionOnEntities
- Parameters
-
y_field_entity_ptr | pointer to second "y_field" field entities |
x_field_entity_ptr | pointer to first "x_field" field entities |
Definition at line 62 of file FieldBlas.hpp.
◆ TwoFieldFunctionOnValues
- Parameters
-
y_field_value_reference | field "y_field" values |
x_field_value | |
Definition at line 53 of file FieldBlas.hpp.
◆ VertexCoordsFunction
◆ FieldBlas()
MoFEM::FieldBlas::FieldBlas |
( |
const MoFEM::Core & |
core | ) |
|
◆ ~FieldBlas()
MoFEM::FieldBlas::~FieldBlas |
( |
| ) |
|
|
default |
◆ fieldAxpy() [1/2]
MoFEMErrorCode MoFEM::FieldBlas::fieldAxpy |
( |
const double |
alpha, |
|
|
const std::string |
field_name_x, |
|
|
const std::string |
field_name_y, |
|
|
bool |
error_if_missing, |
|
|
bool |
create_if_missing |
|
) |
| |
- Deprecated:
- axpy fields
field_y = field_y + alpha*field_x
- Parameters
-
alpha | |
field_name_x | name of field_x |
field_name_y | name of field_y |
error_if_missing | throw error if entity/dof exist in field_x but not on field_y |
create_if_missing | creat dof in field_y from field_x if it is not database |
Definition at line 280 of file FieldBlas.cpp.
285 return fieldAxpy(alpha, field_name_x, field_name_y, error_if_missing);
◆ fieldAxpy() [2/2]
MoFEMErrorCode MoFEM::FieldBlas::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
field_y = field_y + alpha*field_x
- Parameters
-
alpha | |
field_name_x | name of field_x |
field_name_y | name of field_y |
error_if_missing | throw error if entity/dof exist in field_x but not |
ents_ptr | execute lambda function only on entities given by pointer on field_y |
Definition at line 267 of file FieldBlas.cpp.
272 auto axpy = [alpha](
const double fy,
const double fx) {
273 return fy + alpha * fx;
◆ fieldCopy() [1/2]
MoFEMErrorCode MoFEM::FieldBlas::fieldCopy |
( |
const double |
alpha, |
|
|
const std::string |
field_name_x, |
|
|
const std::string |
field_name_y, |
|
|
bool |
error_if_missing, |
|
|
bool |
create_if_missing |
|
) |
| |
◆ fieldCopy() [2/2]
MoFEMErrorCode MoFEM::FieldBlas::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
field_y = alpha*field_x
- Parameters
-
alpha | |
field_name_x | name of field_x |
field_name_y | name of field_y |
error_if_missing | throw error if entity/dof exist in field_x but not |
ents_ptr | execute lambda function only on entities given by pointer on field_y |
create_if_missing | creat dof in field_y from field_x if it is not database |
Definition at line 288 of file FieldBlas.cpp.
293 auto copy = [alpha](
const double fy,
const double fx) {
return alpha * fx; };
◆ fieldLambda()
- Deprecated:
- This is with obsolete last option
- Parameters
-
lambda | |
field_name_x | |
field_name_y | |
error_if_missing | |
create_if_missing | |
- Returns
- DEPRECATED
Definition at line 255 of file FieldBlas.cpp.
260 if (create_if_missing)
262 <<
"Option create_if_missing is set to true but have no effect";
◆ fieldLambdaOnEntities() [1/2]
filed lambda
Do calculation on one field using lambda function with entity field input
auto field_abs = [&](
boost::shared_ptr<FieldEntity> ent_ptr) {
auto field_data = ent_ptr->getEntFieldData();
for (
auto &
v : field_data)
};
"U", block_ents_ptr);
- Parameters
-
function | f(double &x, double) |
field_name | name of field |
ents_ptr | execute lambda function only on entities given by pointer to range |
Definition at line 50 of file FieldBlas.cpp.
58 auto fit = fields_ptr->get<FieldName_mi_tag>().find(
field_name);
59 if (fit == fields_ptr->get<FieldName_mi_tag>().end()) {
61 "field < %s > not found, (top tip: check spelling)",
66 const auto bit_number = (*fit)->getBitNumber();
75 auto eit = field_ents->get<Unique_mi_tag>().lower_bound(lo_uid);
77 auto eit_hi = field_ents->get<Unique_mi_tag>().upper_bound(hi_uid);
79 for (; eit != eit_hi; ++eit) {
87 for (
auto p = ents_ptr->const_pair_begin(); p != ents_ptr->const_pair_end();
89 CHKERR loop_ents(p->first, p->second);
94 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
95 get_id_for_max_type<MBENTITYSET>());
◆ fieldLambdaOnEntities() [2/2]
filed lambda
Do calculation on two fields and save result to field fy input function usees field entities
auto vector_times_scalar_field = [&](boost::shared_ptr<FieldEntity> ent_ptr_y,
boost::shared_ptr<FieldEntity> ent_ptr_x) {
auto x_data = ent_ptr_x->getEntFieldData();
auto y_data = ent_ptr_y->getEntFieldData();
const auto size_x = x_data.size();
const auto size_y = y_data.size();
for (
size_t dd = 0;
dd != size_y; ++
dd)
};
"U", "P" false, block_ents_ptr);
- Parameters
-
function | f(double &x, double) |
field_name_x | name of field_x |
field_name_y | name of field_y |
error_if_missing | throw error if missing entities of field y |
ents_ptr | execute lambda function only on entities given by pointer to range |
Definition at line 158 of file FieldBlas.cpp.
167 auto x_fit = fields_ptr->get<FieldName_mi_tag>().find(field_name_x);
168 if (x_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
170 "x field < %s > not found, (top tip: check spelling)",
171 field_name_x.c_str());
173 auto y_fit = fields_ptr->get<FieldName_mi_tag>().find(field_name_y);
174 if (y_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
176 "y field < %s > not found, (top tip: check spelling)",
177 field_name_y.c_str());
181 const auto x_bit_number = (*x_fit)->getBitNumber();
182 const auto y_bit_number = (*y_fit)->getBitNumber();
183 const auto y_eit_hi = field_ents->get<Unique_mi_tag>().upper_bound(
189 const auto x_lo_uid =
191 const auto x_hi_uid =
195 auto x_eit = field_ents->get<Unique_mi_tag>().lower_bound(x_lo_uid);
197 auto x_eit_hi = field_ents->get<Unique_mi_tag>().upper_bound(x_hi_uid);
199 for (; x_eit != x_eit_hi;) {
202 (*y_fit)->getBitNumber(), (*x_eit)->getEnt());
203 auto y_eit = field_ents->get<Unique_mi_tag>().find(y_lo_uid);
205 if (y_eit == field_ents->end()) {
207 if (error_if_missing) {
209 "Missing entity in y_field.");
217 if (x_eit == x_eit_hi)
219 if (y_eit == y_eit_hi)
221 if ((*y_eit)->getEnt() != (*x_eit)->getEnt())
241 for (
auto p = ents_ptr->const_pair_begin(); p != ents_ptr->const_pair_end();
243 CHKERR loop_ents(p->first, p->second);
248 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
249 get_id_for_max_type<MBENTITYSET>());
◆ fieldLambdaOnValues() [1/2]
filed lambda
Do calculation on one field using lambda function with field value input
auto scale_field = [&](const double val) {
double time = ts_t;
return val * time;
};
"U", false, block_ents_ptr);
- Parameters
-
function | f(const double) |
field_name | name of field |
ents_ptr | execute lambda function only on entities given by pointer to range |
Definition at line 21 of file FieldBlas.cpp.
27 auto x_fit = fields_ptr->get<FieldName_mi_tag>().find(
field_name);
28 if (x_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
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)
◆ fieldLambdaOnValues() [2/2]
filed lambda
Do calculation on two fields and save result to field fy input function usees field values
auto field_axpy = [&](const double val_y,
const double val_x) {
return 2 * val_x + val_y;
};
"U", "P" false, block_ents_ptr);
- Parameters
-
function | f(double &x, double) |
field_name_x | name of field_x |
field_name_y | name of field_y |
error_if_missing | throw error if missing entities of field y |
ents_ptr | execute lambda function only on entities given by pointer to range |
Definition at line 101 of file FieldBlas.cpp.
108 auto x_fit = fields_ptr->get<FieldName_mi_tag>().find(field_name_x);
109 if (x_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
111 "x field < %s > not found, (top tip: check spelling)",
112 field_name_x.c_str());
114 auto y_fit = fields_ptr->get<FieldName_mi_tag>().find(field_name_y);
115 if (y_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
117 "y field < %s > not found, (top tip: check spelling)",
118 field_name_y.c_str());
120 if ((*x_fit)->getSpace() != (*y_fit)->getSpace()) {
122 "space for field < %s > and field <%s> are not compatible for "
124 field_name_x.c_str(), field_name_y.c_str());
126 if ((*x_fit)->getNbOfCoeffs() != (*y_fit)->getNbOfCoeffs()) {
129 "rank for field < %s > and field <%s> are not compatible for fieldblas",
130 field_name_x.c_str(), field_name_y.c_str());
133 auto wrap_lambda_on_enties =
134 [
lambda](boost::shared_ptr<FieldEntity> y_field_entity_ptr,
135 const boost::shared_ptr<FieldEntity> x_field_entity_ptr) {
138 auto x_field_data = x_field_entity_ptr->getEntFieldData();
139 auto y_field_data = y_field_entity_ptr->getEntFieldData();
140 const auto size_x = x_field_data.size();
141 const auto size_y = y_field_data.size();
144 for (;
dd != std::min(size_x, size_y); ++
dd)
145 y_field_data[
dd] =
lambda(y_field_data[
dd], x_field_data[
dd]);
146 for (;
dd < size_y; ++
dd)
147 y_field_data[
dd] = 0;
153 field_name_y, error_if_missing, ents_ptr);
◆ fieldScale()
field scale
- Parameters
-
scale | field by value |
field_name | name of field |
error_if_missing | throw error if missing entities of field y |
ents_ptr | execute lambda function only on entities given by pointer to range |
Definition at line 307 of file FieldBlas.cpp.
312 auto scale = [alpha](
const double v) {
return v * alpha; };
◆ query_interface()
◆ setField() [1/3]
MoFEMErrorCode MoFEM::FieldBlas::setField |
( |
const double |
val, |
|
|
const EntityType |
type, |
|
|
const Range & |
ents, |
|
|
const std::string |
field_name |
|
) |
| |
set field
field_y = val
- Parameters
-
val | |
entity | type |
field_name | |
Definition at line 407 of file FieldBlas.cpp.
420 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(lo_uid);
421 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(hi_uid);
425 for (; dit != hi_dit; dit++) {
426 ent = (*dit)->getEnt();
428 if (ents.find(ent) == ents.end()) {
436 (*dit)->getFieldData() = val;
◆ setField() [2/3]
MoFEMErrorCode MoFEM::FieldBlas::setField |
( |
const double |
val, |
|
|
const EntityType |
type, |
|
|
const std::string |
field_name |
|
) |
| |
- Deprecated:
- scale field
- Parameters
-
val | is a set parameter \field_name is a field name |
Definition at line 386 of file FieldBlas.cpp.
398 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(lo_uid);
399 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(hi_uid);
401 for (; dit != hi_dit; dit++)
402 (*dit)->getFieldData() = val;
◆ setField() [3/3]
set field
field_y = val
- Parameters
-
Definition at line 441 of file FieldBlas.cpp.
447 auto fit = fields_ptr->get<FieldName_mi_tag>().find(
field_name);
448 if (fit == fields_ptr->get<FieldName_mi_tag>().end()) {
450 " field < %s > not found, (top tip: check spelling)",
454 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
456 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
458 for (; dit != hi_dit; dit++) {
459 (*dit)->getFieldData() = val;
◆ setVertexDofs()
Set DOFs on vertices using user function.
Example:
double *y, double *z) {
field_data[0] = (*x);
field_data[1] = (*y);
field_data[2] = (*z);
};
"DISP");
- Note
- Function works both ways, using it coordinates can be set from field.
- Parameters
-
lambda | function evaluating field at points |
field_name | is a field name |
verts | pointer to vertices if null all vertices in the field are evaluated) |
Definition at line 318 of file FieldBlas.cpp.
326 CHKERR m_field.
get_moab().get_entities_by_type(meshset, MBVERTEX, verts,
329 verts = intersect(*sub_verts, verts);
331 struct LambdaMethod : EntityMethod {
334 : EntityMethod(), mField(m_field), verts(verts),
lambda(
lambda),
335 count(0), total(0) {}
342 if (*vit != entPtr->getEnt())
344 "Inconsistent entity %ld != %ld", *vit, entPtr->getEnt());
346 CHKERR mField.get_moab().coords_iterate(vit, verts.end(), xCoord,
347 yCoord, zCoord, count);
348 CHKERR lambda(entPtr->getEntFieldData(), xCoord, yCoord, zCoord);
359 if (total != verts.size())
361 "Inconsistent number of vertices in field meshset and in the "
363 total, verts.size());
◆ cOre
◆ dEbug
bool MoFEM::FieldBlas::dEbug |
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode fieldLambdaOnEntities(OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
MoFEMErrorCode setVertexDofs(VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr)
Set DOFs on vertices using user function.
static UId getLoLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
const MoFEM::Interface & cOre
virtual EntityHandle get_field_meshset(const std::string name) const =0
get field meshset
Deprecated interface functions.
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
#define CHKERR
Inline error check.
virtual moab::Interface & get_moab()=0
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
EntityHandle get_id_for_max_type()
virtual const FieldEntity_multiIndex * get_field_ents() const =0
Get the field ents object.
VectorShallowArrayAdaptor< double > VectorAdaptor
constexpr auto field_name
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
virtual const Field_multiIndex * get_fields() const =0
Get the fields object.
const double v
phase velocity of light in medium (cm/ns)
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)
#define MOFEM_LOG(channel, severity)
Log.
FieldBlas(const MoFEM::Core &core)
@ MOFEM_DATA_INCONSISTENCY
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
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static UId getHiLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
boost::function< MoFEMErrorCode(VectorAdaptor &&field_data, double *xcoord, double *ycoord, double *zcoord)> VertexCoordsFunction
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
EntityHandle get_id_for_min_type()
#define MoFEMFunctionReturn(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 ...