v0.14.0
Loading...
Searching...
No Matches
FieldBlas.cpp
Go to the documentation of this file.
1/** \file FieldBlas.cpp
2 * \brief Managing complexities for problem
3 * \ingroup mofem_section_manager
4 */
5
6
7namespace MoFEM {
8
10FieldBlas::query_interface(boost::typeindex::type_index type_index,
11 UnknownInterface **iface) const {
13 *iface = const_cast<FieldBlas *>(this);
15}
16
18 : cOre(const_cast<MoFEM::Core &>(core)), dEbug(false) {}
19
22 const std::string field_name, Range *ents_ptr) {
23 const MoFEM::Interface &m_field = cOre;
24 auto fields_ptr = m_field.get_fields();
26
27 auto x_fit = fields_ptr->get<FieldName_mi_tag>().find(field_name);
28 if (x_fit == fields_ptr->get<FieldName_mi_tag>().end()) {
29 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
30 "field < %s > not found, (top tip: check spelling)",
31 field_name.c_str());
32 }
33
34 auto wrap_lambda_on_enties =
35 [lambda](boost::shared_ptr<FieldEntity> field_entity_ptr) {
37
38 auto field_data = field_entity_ptr->getEntFieldData();
39 for (auto &v : field_data)
40 v = lambda(v);
41
43 };
44
45 CHKERR fieldLambdaOnEntities(wrap_lambda_on_enties, field_name, ents_ptr);
47};
48
51 const std::string field_name,
52 Range *ents_ptr) {
53 const MoFEM::Interface &m_field = cOre;
54 auto fields_ptr = m_field.get_fields();
55 auto field_ents = m_field.get_field_ents();
57
58 auto fit = fields_ptr->get<FieldName_mi_tag>().find(field_name);
59 if (fit == fields_ptr->get<FieldName_mi_tag>().end()) {
60 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
61 "field < %s > not found, (top tip: check spelling)",
62 field_name.c_str());
63 }
64
65 // End of iterator for y_field
66 const auto bit_number = (*fit)->getBitNumber();
67
68 auto loop_ents = [&](const EntityHandle f, const EntityHandle s) {
70
71 const auto lo_uid = FieldEntity::getLoLocalEntityBitNumber(bit_number, f);
72 const auto hi_uid = FieldEntity::getHiLocalEntityBitNumber(bit_number, s);
73
74 // Start of iterator for x_field
75 auto eit = field_ents->get<Unique_mi_tag>().lower_bound(lo_uid);
76 // End of iterator for x_field
77 auto eit_hi = field_ents->get<Unique_mi_tag>().upper_bound(hi_uid);
78
79 for (; eit != eit_hi; ++eit) {
80 CHKERR lambda(*eit);
81 }
82
84 };
85
86 if (ents_ptr) {
87 for (auto p = ents_ptr->const_pair_begin(); p != ents_ptr->const_pair_end();
88 ++p) {
89 CHKERR loop_ents(p->first, p->second);
90 }
91 } else {
92 // we are looping from the very first possible entity handle (MBVERTEX) to
93 // the very last (MBENTITYSET)
94 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
95 get_id_for_max_type<MBENTITYSET>());
96 }
97
99}
100
102 FieldBlas::TwoFieldFunctionOnValues lambda, const std::string field_name_x,
103 const std::string field_name_y, bool error_if_missing, Range *ents_ptr) {
104 const MoFEM::Interface &m_field = cOre;
105 auto fields_ptr = m_field.get_fields();
107
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()) {
110 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
111 "x field < %s > not found, (top tip: check spelling)",
112 field_name_x.c_str());
113 }
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()) {
116 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
117 "y field < %s > not found, (top tip: check spelling)",
118 field_name_y.c_str());
119 }
120 if ((*x_fit)->getSpace() != (*y_fit)->getSpace()) {
121 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
122 "space for field < %s > and field <%s> are not compatible for "
123 "fieldblas",
124 field_name_x.c_str(), field_name_y.c_str());
125 }
126 if ((*x_fit)->getNbOfCoeffs() != (*y_fit)->getNbOfCoeffs()) {
127 SETERRQ2(
128 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
129 "rank for field < %s > and field <%s> are not compatible for fieldblas",
130 field_name_x.c_str(), field_name_y.c_str());
131 }
132
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) {
137
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();
142
143 size_t dd = 0;
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;
148
150 };
151
152 CHKERR fieldLambdaOnEntities(wrap_lambda_on_enties, field_name_x,
153 field_name_y, error_if_missing, ents_ptr);
155};
156
159 const std::string field_name_x,
160 const std::string field_name_y,
161 bool error_if_missing, Range *ents_ptr) {
162 const MoFEM::Interface &m_field = cOre;
163 auto fields_ptr = m_field.get_fields();
164 auto field_ents = m_field.get_field_ents();
166
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()) {
169 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
170 "x field < %s > not found, (top tip: check spelling)",
171 field_name_x.c_str());
172 }
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()) {
175 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
176 "y field < %s > not found, (top tip: check spelling)",
177 field_name_y.c_str());
178 }
179
180 // End of iterator for y_field
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(
184 FieldEntity::getHiBitNumberUId(y_bit_number));
185
186 auto loop_ents = [&](const EntityHandle f, const EntityHandle s) {
188
189 const auto x_lo_uid =
191 const auto x_hi_uid =
193
194 // Start of iterator for x_field
195 auto x_eit = field_ents->get<Unique_mi_tag>().lower_bound(x_lo_uid);
196 // End of iterator for x_field
197 auto x_eit_hi = field_ents->get<Unique_mi_tag>().upper_bound(x_hi_uid);
198
199 for (; x_eit != x_eit_hi;) {
200
201 const auto y_lo_uid = FieldEntity::getLocalUniqueIdCalculate(
202 (*y_fit)->getBitNumber(), (*x_eit)->getEnt());
203 auto y_eit = field_ents->get<Unique_mi_tag>().find(y_lo_uid);
204
205 if (y_eit == field_ents->end()) {
206
207 if (error_if_missing) {
208 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
209 "Missing entity in y_field.");
210 }
211
212 ++x_eit;
213
214 } else {
215
216 auto check = [&]() {
217 if (x_eit == x_eit_hi)
218 return false;
219 if (y_eit == y_eit_hi)
220 return false;
221 if ((*y_eit)->getEnt() != (*x_eit)->getEnt())
222 return false;
223 return true;
224 };
225
226 do {
227
228 CHKERR lambda(*y_eit, *x_eit);
229
230 ++x_eit;
231 ++y_eit;
232
233 } while (check());
234 }
235 }
236
238 };
239
240 if (ents_ptr) {
241 for (auto p = ents_ptr->const_pair_begin(); p != ents_ptr->const_pair_end();
242 ++p) {
243 CHKERR loop_ents(p->first, p->second);
244 }
245 } else {
246 // we are looping from the very first possible entity handle (MBVERTEX) to
247 // the very last (MBENTITYSET)
248 CHKERR loop_ents(get_id_for_min_type<MBVERTEX>(),
249 get_id_for_max_type<MBENTITYSET>());
250 }
251
253}
254
256 const std::string field_name_x,
257 const std::string field_name_y,
258 bool error_if_missing,
259 bool create_if_missing) {
260 if (create_if_missing)
261 MOFEM_LOG("SELF", Sev::noisy)
262 << "Option create_if_missing is set to true but have no effect";
263 return fieldLambdaOnValues(lambda, field_name_x, field_name_y,
264 error_if_missing);
265}
266
268 const std::string field_name_x,
269 const std::string field_name_y,
270 bool error_if_missing, Range *ents_ptr) {
272 auto axpy = [alpha](const double fy, const double fx) {
273 return fy + alpha * fx;
274 };
275 CHKERR fieldLambdaOnValues(axpy, field_name_x, field_name_y, error_if_missing,
276 ents_ptr);
278}
279
281 const std::string field_name_x,
282 const std::string field_name_y,
283 bool error_if_missing,
284 bool create_if_missing) {
285 return fieldAxpy(alpha, field_name_x, field_name_y, error_if_missing);
286}
287
289 const std::string field_name_x,
290 const std::string field_name_y,
291 bool error_if_missing, Range *ents_ptr) {
293 auto copy = [alpha](const double fy, const double fx) { return alpha * fx; };
294 CHKERR fieldLambdaOnValues(copy, field_name_x, field_name_y, error_if_missing,
295 ents_ptr);
297}
298
300 const std::string field_name_x,
301 const std::string field_name_y,
302 bool error_if_missing,
303 bool create_if_missing) {
304 return fieldCopy(alpha, field_name_x, field_name_y, error_if_missing);
305}
306
308 const std::string field_name,
309 Range *ents_ptr) {
311
312 auto scale = [alpha](const double v) { return v * alpha; };
314
316}
317
319 const std::string field_name,
320 Range *sub_verts) {
321 const MoFEM::Interface &m_field = cOre;
323
324 EntityHandle meshset = m_field.get_field_meshset(field_name);
325 Range verts;
326 CHKERR m_field.get_moab().get_entities_by_type(meshset, MBVERTEX, verts,
327 true);
328 if (sub_verts)
329 verts = intersect(*sub_verts, verts);
330
331 struct LambdaMethod : EntityMethod {
332 LambdaMethod(MoFEM::Interface &m_field, Range &verts,
334 : EntityMethod(), mField(m_field), verts(verts), lambda(lambda),
335 count(0), total(0) {}
336 MoFEMErrorCode preProcess() {
337 vit = verts.begin();
338 return 0;
339 }
340 MoFEMErrorCode operator()() {
342 if (*vit != entPtr->getEnt())
343 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
344 "Inconsistent entity %ld != %ld", *vit, entPtr->getEnt());
345 if (!count)
346 CHKERR mField.get_moab().coords_iterate(vit, verts.end(), xCoord,
347 yCoord, zCoord, count);
348 CHKERR lambda(entPtr->getEntFieldData(), xCoord, yCoord, zCoord);
349 ++xCoord;
350 ++yCoord;
351 ++zCoord;
352 ++vit;
353 ++total;
354 --count;
356 }
357 MoFEMErrorCode postProcess() {
359 if (total != verts.size())
360 SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
361 "Inconsistent number of vertices in field meshset and in the "
362 "field %d != %d",
363 total, verts.size());
365 }
366
367 private:
368 MoFEM::Interface &mField;
369 Range &verts;
371 int count;
372 int total;
373 Range::iterator vit;
374 double *xCoord;
375 double *yCoord;
376 double *zCoord;
377 };
378
379 LambdaMethod lambda_method(const_cast<MoFEM::Interface &>(m_field), verts,
380 lambda);
381 CHKERR const_cast<MoFEM::Interface &>(m_field).loop_entities(
382 field_name, lambda_method, &verts, QUIET);
384}
385
386MoFEMErrorCode FieldBlas::setField(const double val, const EntityType type,
387 const std::string field_name) {
388 const MoFEM::Interface &m_field = cOre;
389 auto dofs_ptr = m_field.get_dofs();
391
392 const auto bit_number = m_field.get_field_bit_number(field_name);
393 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
394 bit_number, get_id_for_min_type(type));
395 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
396 bit_number, get_id_for_max_type(type));
397
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);
400
401 for (; dit != hi_dit; dit++)
402 (*dit)->getFieldData() = val;
403
405}
406
407MoFEMErrorCode FieldBlas::setField(const double val, const EntityType type,
408 const Range &ents,
409 const std::string field_name) {
410 const MoFEM::Interface &m_field = cOre;
411 auto dofs_ptr = m_field.get_dofs();
413
414 const auto bit_number = m_field.get_field_bit_number(field_name);
415 const auto lo_uid = FieldEntity::getLocalUniqueIdCalculate(
416 bit_number, get_id_for_min_type(type));
417 const auto hi_uid = FieldEntity::getLocalUniqueIdCalculate(
418 bit_number, get_id_for_max_type(type));
419
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);
422
423 EntityHandle ent, last = 0;
424 bool cont = true;
425 for (; dit != hi_dit; dit++) {
426 ent = (*dit)->getEnt();
427 if (ent != last) {
428 if (ents.find(ent) == ents.end()) {
429 cont = true;
430 } else {
431 cont = false;
432 }
433 last = ent;
434 }
435 if (!cont)
436 (*dit)->getFieldData() = val;
437 }
439}
440
442 const std::string field_name) {
443 const MoFEM::Interface &m_field = cOre;
444 auto fields_ptr = m_field.get_fields();
445 auto dofs_ptr = m_field.get_dofs();
447 auto fit = fields_ptr->get<FieldName_mi_tag>().find(field_name);
448 if (fit == fields_ptr->get<FieldName_mi_tag>().end()) {
449 SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
450 " field < %s > not found, (top tip: check spelling)",
451 field_name.c_str());
452 }
453
454 auto dit = dofs_ptr->get<Unique_mi_tag>().lower_bound(
455 FieldEntity::getLoBitNumberUId((*fit)->getBitNumber()));
456 auto hi_dit = dofs_ptr->get<Unique_mi_tag>().upper_bound(
457 FieldEntity::getHiBitNumberUId((*fit)->getBitNumber()));
458 for (; dit != hi_dit; dit++) {
459 (*dit)->getFieldData() = val;
460 }
462}
463
464} // namespace MoFEM
static Index< 'p', 3 > p
@ QUIET
Definition: definitions.h:208
#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
#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
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.
Definition: LogManager.hpp:308
static double lambda
const double v
phase velocity of light in medium (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
EntityHandle get_id_for_max_type()
EntityHandle get_id_for_min_type()
double scale
Definition: plastic.cpp:170
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
Core (interface) class.
Definition: Core.hpp:82
Deprecated interface functions.
Data structure to exchange data between mofem and User Loop Methods on entities.
Basic algebra on fields.
Definition: FieldBlas.hpp:21
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)
Definition: FieldBlas.cpp:255
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
Definition: FieldBlas.cpp:267
boost::function< MoFEMErrorCode(boost::shared_ptr< FieldEntity > field_entity_ptr)> OneFieldFunctionOnEntities
Definition: FieldBlas.hpp:45
MoFEMErrorCode fieldLambdaOnEntities(OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
Definition: FieldBlas.cpp:50
boost::function< double(const double y_field_value_reference, const double x_field_value)> TwoFieldFunctionOnValues
Definition: FieldBlas.hpp:53
MoFEMErrorCode setVertexDofs(VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr)
Set DOFs on vertices using user function.
Definition: FieldBlas.cpp:318
boost::function< MoFEMErrorCode(VectorAdaptor &&field_data, double *xcoord, double *ycoord, double *zcoord)> VertexCoordsFunction
Definition: FieldBlas.hpp:282
MoFEMErrorCode setField(const double val, const EntityType type, const std::string field_name)
Definition: FieldBlas.cpp:386
boost::function< double(const double field_val)> OneFieldFunctionOnValues
function to set a field value
Definition: FieldBlas.hpp:38
MoFEMErrorCode fieldScale(const double alpha, const std::string field_name, Range *ents_ptr=nullptr)
field scale
Definition: FieldBlas.cpp:307
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: FieldBlas.cpp:10
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
Definition: FieldBlas.cpp:21
FieldBlas(const MoFEM::Core &core)
Definition: FieldBlas.cpp:17
const MoFEM::Interface & cOre
Definition: FieldBlas.hpp:26
boost::function< MoFEMErrorCode(boost::shared_ptr< FieldEntity > y_field_entity_ptr, const boost::shared_ptr< FieldEntity > x_field_entity_ptr)> TwoFieldFunctionOnEntities
Definition: FieldBlas.hpp:62
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
Definition: FieldBlas.cpp:288
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