v0.14.0
FieldBlas.hpp
Go to the documentation of this file.
1 /** \file FieldBlas.hpp
2  * \brief Field basic algebra
3  * \ingroup mofem_is_managers
4  *
5  * Managing problems, build and partitioning.
6  *
7  */
8 
9 #ifndef __FIELD_BLAS_HPP__
10 #define __FIELD_BLAS_HPP__
11 
12 #include "UnknownInterface.hpp"
13 
14 namespace MoFEM {
15 
16 /**
17  * \brief Basic algebra on fields
18  * \ingroup mofem_field_algebra
19  *
20  */
21 struct FieldBlas : public UnknownInterface {
22 
23  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
24  UnknownInterface **iface) const;
25 
27  bool dEbug;
28 
29  FieldBlas(const MoFEM::Core &core);
30 
31  ~FieldBlas() = default;
32 
33  /**
34  * @brief function to set a field value
35  *
36  */
38  boost::function<double(const double field_val)>;
39 
40  /**
41  * @brief
42  *
43  */
44  using OneFieldFunctionOnEntities = boost::function<MoFEMErrorCode(
45  boost::shared_ptr<FieldEntity> field_entity_ptr)>;
46 
47  /**
48  * @brief
49  * \param y_field_value_reference field "y_field" values
50  * \param x_field_value
51  */
52  using TwoFieldFunctionOnValues = boost::function<double(
53  const double y_field_value_reference, const double x_field_value)>;
54 
55  /**
56  * @brief
57  * \param y_field_entity_ptr pointer to second "y_field" field entities
58  * \param x_field_entity_ptr pointer to first "x_field" field entities
59  */
60  using TwoFieldFunctionOnEntities = boost::function<MoFEMErrorCode(
61  boost::shared_ptr<FieldEntity> y_field_entity_ptr,
62  const boost::shared_ptr<FieldEntity> x_field_entity_ptr)>;
63 
64  /** \brief filed lambda
65  * \ingroup mofem_field_algebra
66  *
67  * Do calculation on one field using lambda function with entity field input
68  *
69  * \code
70  auto field_abs = [&](
71  boost::shared_ptr<FieldEntity> ent_ptr) {
72  MoFEMFunctionBeginHot;
73  auto field_data = ent_ptr->getEntFieldData();
74  for (auto &v : field_data)
75  v = std::abs(v);
76  MoFEMFunctionReturnHot(0);
77  };
78 
79  CHKERR
80  m_field.getInterface<FieldBlas>()->fieldLambdaOnEntities(field_abs,
81  "U", block_ents_ptr);
82  * \endcode
83  *
84  * \param function f(double &x, double)
85  * \param field_name name of field
86  * \param ents_ptr execute lambda function only on entities given by pointer
87  to range
88  *
89  */
91  const std::string field_name,
92  Range *ents_ptr = nullptr);
93 
94  /** \brief filed lambda
95  * \ingroup mofem_field_algebra
96  *
97  * Do calculation on one field using lambda function with field value input
98  *
99  * \code
100  auto scale_field = [&](const double val) {
101  double time = ts_t; // global variable
102  return val * time;
103  };
104  CHKERR
105  m_field.getInterface<FieldBlas>()->fieldLambdaOnValues(scale_field,
106  "U", false, block_ents_ptr);
107  * \endcode
108  *
109  * \param function f(const double)
110  * \param field_name name of field
111  * \param ents_ptr execute lambda function only on entities given by
112  pointer to range
113  *
114  */
116  const std::string field_name,
117  Range *ents_ptr = nullptr);
118 
119  /** \brief filed lambda
120  * \ingroup mofem_field_algebra
121  *
122  * Do calculation on two fields and save result to field fy
123  * input function usees field values
124  * \code
125  auto field_axpy = [&](const double val_y,
126  const double val_x) {
127  // y = 2 * x + y
128  return 2 * val_x + val_y;
129  };
130  CHKERR
131  m_field.getInterface<FieldBlas>()->fieldLambdaOnValues(field_axpy,
132  "U", "P" false, block_ents_ptr);
133  * \endcode
134  *
135  * \param function f(double &x, double)
136  * \param field_name_x name of field_x
137  * \param field_name_y name of field_y
138  * \param error_if_missing throw error if missing entities of field y
139  * \param ents_ptr execute lambda function only on entities given by pointer
140  to range
141  *
142  */
144  const std::string field_name_x,
145  const std::string field_name_y,
146  bool error_if_missing = false,
147  Range *ents_ptr = nullptr);
148 
149  /** \brief filed lambda
150  * \ingroup mofem_field_algebra
151  *
152  * Do calculation on two fields and save result to field fy
153  * input function usees field entities
154  *
155  * \code
156  auto vector_times_scalar_field = [&](boost::shared_ptr<FieldEntity> ent_ptr_y,
157  boost::shared_ptr<FieldEntity> ent_ptr_x) {
158  MoFEMFunctionBeginHot;
159  auto x_data = ent_ptr_x->getEntFieldData();
160  auto y_data = ent_ptr_y->getEntFieldData();
161  const auto size_x = x_data.size(); // scalar
162  const auto size_y = y_data.size(); // vector
163 
164  for (size_t dd = 0; dd != size_y; ++dd)
165  y_data[dd] *= x_data[0];
166 
167  MoFEMFunctionReturnHot(0);
168  };
169 
170  CHKERR
171  m_field.getInterface<FieldBlas>()->fieldLambdaOnValues(vector_times_scalar_field,
172  "U", "P" false, block_ents_ptr);
173 
174  * \endcode
175  *
176  * \param function f(double &x, double)
177  * \param field_name_x name of field_x
178  * \param field_name_y name of field_y
179  * \param error_if_missing throw error if missing entities of field y
180  * \param ents_ptr execute lambda function only on entities given by pointer
181  to range
182  *
183  */
185  const std::string field_name_x,
186  const std::string field_name_y,
187  bool error_if_missing = false,
188  Range *ents_ptr = nullptr);
189 
190  /**
191  * @deprecated This is with obsolete last option
192  *
193  * @param lambda
194  * @param field_name_x
195  * @param field_name_y
196  * @param error_if_missing
197  * @param create_if_missing
198  * @return DEPRECATED
199  */
201  const std::string field_name_x,
202  const std::string field_name_y,
203  bool error_if_missing,
204  bool create_if_missing);
205 
206  /** \brief axpy fields
207  * \ingroup mofem_field_algebra
208  *
209  * field_y = field_y + alpha*field_x
210  *
211  * \param alpha
212  * \param field_name_x name of field_x
213  * \param field_name_y name of field_y
214  * \param error_if_missing throw error if entity/dof exist in field_x but not
215  * \param ents_ptr execute lambda function only on entities given by pointer
216  * on field_y
217  *
218  */
219  MoFEMErrorCode fieldAxpy(const double alpha, const std::string field_name_x,
220  const std::string field_name_y,
221  bool error_if_missing = false,
222  Range *ents_ptr = nullptr);
223 
224  /** \deprecated axpy fields
225  * \ingroup mofem_field_algebra
226  *
227  * field_y = field_y + alpha*field_x
228  *
229  * \param alpha
230  * \param field_name_x name of field_x
231  * \param field_name_y name of field_y
232  * \param error_if_missing throw error if entity/dof exist in field_x but not
233  * on field_y \param create_if_missing creat dof in field_y from field_x if it
234  * is not database
235  *
236  */
237  DEPRECATED MoFEMErrorCode fieldAxpy(const double alpha,
238  const std::string field_name_x,
239  const std::string field_name_y,
240  bool error_if_missing,
241  bool create_if_missing);
242 
243  /** \brief copy and scale fields
244  * \ingroup mofem_field_algebra
245  *
246  * field_y = alpha*field_x
247  *
248  * \param alpha
249  * \param field_name_x name of field_x
250  * \param field_name_y name of field_y
251  * \param error_if_missing throw error if entity/dof exist in field_x but not
252  * \param ents_ptr execute lambda function only on entities given by pointer
253  * on field_y \param create_if_missing creat dof in field_y from field_x if it
254  * is not database
255  *
256  */
257  MoFEMErrorCode fieldCopy(const double alpha, const std::string field_name_x,
258  const std::string field_name_y,
259  bool error_if_missing = false,
260  Range *ents_ptr = nullptr);
261 
262  DEPRECATED MoFEMErrorCode fieldCopy(const double alpha,
263  const std::string field_name_x,
264  const std::string field_name_y,
265  bool error_if_missing,
266  bool create_if_missing);
267 
268  /** \brief field scale
269  * \ingroup mofem_field_algebra
270  *
271  * \param scale field by value
272  * \param field_name name of field
273  * \param error_if_missing throw error if missing entities of field y
274  * \param ents_ptr execute lambda function only on entities given by
275  pointer to range
276  */
277  MoFEMErrorCode fieldScale(const double alpha, const std::string field_name,
278  Range *ents_ptr = nullptr);
279 
280  using VertexCoordsFunction =
281  boost::function<MoFEMErrorCode(VectorAdaptor &&field_data, double *xcoord,
282  double *ycoord, double *zcoord)>;
283 
284  /** \brief Set DOFs on vertices using user function
285  * \ingroup mofem_field_algebra
286  *
287  * Example:
288  *
289  * \code
290  * auto do_something = [&](VectorAdaptor &field_data, double *x,
291  * double *y, double *z) {
292  * MoFEMFunctionBegin;
293  * field_data[0] = (*x);
294  * field_data[1] = (*y);
295  * field_data[2] = (*z);
296  * MoFEMFunctionReturn(0);
297  * };
298  * CHKERR m_field.getInterface<FieldBlas>()->setVertexDofs(set_distance,
299  * "DISP"); \endcode
300  *
301  * \note Function works both ways, using it coordinates can be set from
302  * field.
303  *
304  * \param lambda function evaluating field at points
305  * \param field_name is a field name
306  * \param verts pointer to vertices if null all vertices in the field are
307  * evaluated)
308  *
309  */
311  const std::string field_name,
312  Range *verts = nullptr);
313 
314  /** \deprecated scale field
315  * \ingroup mofem_field_algebra
316  *
317  * \param val is a set parameter
318  * \field_name is a field name
319  *
320  */
321  MoFEMErrorCode setField(const double val, const EntityType type,
322  const std::string field_name);
323 
324  /** \brief set field
325  * \ingroup mofem_field_algebra
326  *
327  * field_y = val
328  *
329  * \param val
330  * \param entity type
331  * \param field_name
332  *
333  */
334  MoFEMErrorCode setField(const double val, const EntityType type,
335  const Range &ents, const std::string field_name);
336 
337  /** \brief set field
338  * \ingroup mofem_field_algebra
339  *
340  * field_y = val
341  *
342  * \param val
343  * \param field_name
344  *
345  */
346  MoFEMErrorCode setField(const double val, const std::string field_name);
347 
348 };
349 
350 } // namespace MoFEM
351 
352 #endif // __FIELD_BLAS_HPP__
MoFEM::FieldBlas::fieldLambdaOnEntities
MoFEMErrorCode fieldLambdaOnEntities(OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
Definition: FieldBlas.cpp:50
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::FieldBlas::setVertexDofs
MoFEMErrorCode setVertexDofs(VertexCoordsFunction lambda, const std::string field_name, Range *verts=nullptr)
Set DOFs on vertices using user function.
Definition: FieldBlas.cpp:318
DEPRECATED
#define DEPRECATED
Definition: definitions.h:17
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::FieldBlas::cOre
const MoFEM::Interface & cOre
Definition: FieldBlas.hpp:26
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::FieldBlas::fieldLambda
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
MoFEM::FieldBlas::dEbug
bool dEbug
Definition: FieldBlas.hpp:27
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::FieldBlas::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: FieldBlas.cpp:10
MoFEM::FieldBlas::fieldLambdaOnValues
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
filed lambda
Definition: FieldBlas.cpp:21
MoFEM::FieldBlas::TwoFieldFunctionOnEntities
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
double
convert.type
type
Definition: convert.py:64
MoFEM::FieldBlas::TwoFieldFunctionOnValues
boost::function< double(const double y_field_value_reference, const double x_field_value)> TwoFieldFunctionOnValues
Definition: FieldBlas.hpp:53
UnknownInterface.hpp
MoFEM interface.
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::FieldBlas::fieldCopy
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
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
Range
MoFEM::FieldBlas::~FieldBlas
~FieldBlas()=default
MoFEM::FieldBlas::FieldBlas
FieldBlas(const MoFEM::Core &core)
Definition: FieldBlas.cpp:17
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
MoFEM::FieldBlas::OneFieldFunctionOnEntities
boost::function< MoFEMErrorCode(boost::shared_ptr< FieldEntity > field_entity_ptr)> OneFieldFunctionOnEntities
Definition: FieldBlas.hpp:45
MoFEM::FieldBlas::setField
MoFEMErrorCode setField(const double val, const EntityType type, const std::string field_name)
Definition: FieldBlas.cpp:386
MoFEM::FieldBlas::fieldAxpy
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
MoFEM::FieldBlas::VertexCoordsFunction
boost::function< MoFEMErrorCode(VectorAdaptor &&field_data, double *xcoord, double *ycoord, double *zcoord)> VertexCoordsFunction
Definition: FieldBlas.hpp:282
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
MoFEM::FieldBlas::OneFieldFunctionOnValues
boost::function< double(const double field_val)> OneFieldFunctionOnValues
function to set a field value
Definition: FieldBlas.hpp:38
MoFEM::FieldBlas::fieldScale
MoFEMErrorCode fieldScale(const double alpha, const std::string field_name, Range *ents_ptr=nullptr)
field scale
Definition: FieldBlas.cpp:307