v0.14.0
field_blas_axpy.cpp
Go to the documentation of this file.
1 /** \file field_blas.cpp
2  * \example field_blas.cpp
3  * \brief test field blas interface
4  *
5  * \ingroup mofem_field_algebra
6  */
7 
8 
9 
10 #include <MoFEM.hpp>
11 
12 using namespace MoFEM;
13 
14 static char help[] = "...\n\n";
15 
16 int main(int argc, char *argv[]) {
17 
18  MoFEM::Core::Initialize(&argc, &argv, PETSC_NULL, help);
19 
20  try {
21 
22  moab::Core mb_instance;
23  moab::Interface &moab = mb_instance;
24  int rank;
25  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
26 
27  // Reade parameters from line command
28  PetscBool flg = PETSC_TRUE;
29  char mesh_file_name[255];
30 #if PETSC_VERSION_GE(3, 6, 4)
31  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
32  255, &flg);
33 #else
34  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
35  mesh_file_name, 255, &flg);
36 #endif
37  if (flg != PETSC_TRUE) {
38  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
39  }
40  int order = 1;
41 #if PETSC_VERSION_GE(3, 6, 4)
42  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-my_order", &order, PETSC_NULL);
43 #else
44  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
45  PETSC_NULL);
46 #endif
47 
48  // Read mesh to MOAB
49  const char *option;
50  option = "";
51  CHKERR moab.load_file(mesh_file_name, 0, option);
52 
53  // Create MoFEM database
54  // Note: Is MoFEM::CoreTmp<1> for testing purposes only
55  MoFEM::CoreTmp<-1> core(moab);
56  MoFEM::Interface &m_field = core;
57 
58  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
59  0, 3, BitRefLevel().set(0));
60 
61  // Create fields, add entities to field and set approximation order
62  CHKERR m_field.add_field("FIELD_A", H1, AINSWORTH_LEGENDRE_BASE, 3,
63  MB_TAG_DENSE);
64  CHKERR m_field.add_field("FIELD_B", H1, AINSWORTH_LEGENDRE_BASE, 3,
65  MB_TAG_DENSE);
66  CHKERR m_field.add_field("FIELD_C", H1, AINSWORTH_LEGENDRE_BASE, 1,
67  MB_TAG_DENSE);
68 
69  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_A");
70  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_B");
71  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_C");
72 
73  CHKERR m_field.set_field_order(0, MBTET, "FIELD_A", order + 1);
74  CHKERR m_field.set_field_order(0, MBTRI, "FIELD_A", order + 1);
75  CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_A", order + 1);
76  CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_A", 1);
77 
78  CHKERR m_field.set_field_order(0, MBTET, "FIELD_B", order);
79  CHKERR m_field.set_field_order(0, MBTRI, "FIELD_B", order);
80  CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_B", order);
81  CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_B", 1);
82 
83  CHKERR m_field.set_field_order(0, MBTET, "FIELD_C", order);
84  CHKERR m_field.set_field_order(0, MBTRI, "FIELD_C", order);
85  CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_C", order);
86  CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_C", 1);
87 
88  // build field
89  CHKERR m_field.build_fields();
90 
91  // get access to field agebra interface
92  auto fb = m_field.getInterface<FieldBlas>();
93 
94  auto check_axpy = [&]() {
96  // set value to field
97  CHKERR fb->setField(+1, MBVERTEX, "FIELD_A");
98  CHKERR fb->setField(-2, MBVERTEX, "FIELD_B");
99 
100  // FIELD_A = FIELD_A + 0.5 * FIELD_B
101  CHKERR fb->fieldAxpy(+0.5, "FIELD_B", "FIELD_A");
102  // FIELD_B *= -0.5
103  CHKERR fb->fieldScale(-0.5, "FIELD_B");
104 
105  auto dofs_ptr = m_field.get_dofs();
106  for (auto dit : *dofs_ptr) {
107 
108  auto check = [&](const std::string name, const double expected) {
110  if (dit->getName() == name) {
111  cout << name << " " << dit->getFieldData() << " " << expected
112  << endl;
113  if (dit->getFieldData() != expected)
114  SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
115  "Wrong DOF value 0 != %4.3e for %s", dit->getFieldData(),
116  boost::lexical_cast<std::string>(*dit).c_str());
117  }
118 
120  };
121 
122  CHKERR check("FIELD_A", 0);
123  CHKERR check("FIELD_B", 1);
124  }
126  };
127 
128  auto check_set_vertex = [&]() {
130 
131  auto set_distance = [&](VectorAdaptor &&field_data, double *xcoord,
132  double *ycoord, double *zcoord) {
136  xcoord, ycoord, zcoord);
137  field_data[0] = sqrt(t_coord(i) * t_coord(i));
139  };
140 
141  CHKERR fb->setVertexDofs(set_distance, "FIELD_A");
142 
143  // Test set veryex
144  struct TestMethod : EntityMethod {
145  TestMethod(moab::Interface &moab) : EntityMethod(), moab(moab) {}
146 
147  MoFEMErrorCode preProcess() { return 0; }
148  MoFEMErrorCode operator()() {
150  if (entPtr->getEntType() == MBVERTEX) {
151  EntityHandle v = entPtr->getEnt();
152  VectorDouble3 coords(3);
153  CHKERR moab.get_coords(&v, 1, &coords[0]);
154  if (std::abs(entPtr->getEntFieldData()[0] - norm_2(coords)) >
155  std::numeric_limits<double>::epsilon())
156  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
157  "Wrong vals %3.4e != %3.4e",
158  entPtr->getEntFieldData()[0], norm_2(coords));
159  }
161  }
162  MoFEMErrorCode postProcess() { return 0; }
163 
164  private:
165  moab::Interface &moab;
166  };
167 
168  TestMethod method(moab);
169  // checking if all is ok
170  CHKERR m_field.loop_entities("FIELD_A", method);
171 
173  };
174 
175  auto check_lambda = [&]() {
177 
178  Range meshset_ents;
180  m_field, SIDESET | PRESSURESET, it)) {
181  if (it->getMeshsetId() != 1)
182  continue;
183  Range ents, nodes;
184  CHKERR it->getMeshsetIdEntitiesByType(m_field.get_moab(), MBTRI, ents,
185  true);
186  CHKERR m_field.get_moab().get_connectivity(ents, nodes, true);
187  meshset_ents.merge(nodes);
188  }
189  // set value to field
190  CHKERR fb->setField(+1, MBVERTEX, "FIELD_A");
191  CHKERR fb->setField(-1, MBVERTEX, "FIELD_B");
192  CHKERR fb->setField(1.5, MBVERTEX, "FIELD_C");
193 
194  auto field_axpy = [&](const double val_y, const double val_x) {
195  return 0.5 * val_x + val_y;
196  };
197  auto field_scale = [&](const double val) { return -0.5 * val; };
198 
199  // note that y_ent is first
200  // FIXME: this can be confusing?
201  auto vector_times_scalar_field =
202  [&](boost::shared_ptr<FieldEntity> ent_ptr_y,
203  boost::shared_ptr<FieldEntity> ent_ptr_x) {
205  auto x_data = ent_ptr_x->getEntFieldData();
206  auto y_data = ent_ptr_y->getEntFieldData();
207  const auto size_x = x_data.size(); // scalar
208  const auto size_y = y_data.size(); // vector
209 
210  for (size_t dd = 0; dd != size_y; ++dd)
211  y_data[dd] *= x_data[0];
212 
213  // y_data *= x_data[0]; //would work as well
214 
216  };
217 
218  Range ents; // TODO: create test for subentities
219  // FIELD_A = FIELD_A + 0.5 * FIELD_B
220  CHKERR fb->fieldLambdaOnValues(field_axpy, "FIELD_B", "FIELD_A", true);
221  // CHKERR fb->fieldAxpy(+0.5, "FIELD_B", "FIELD_A");
222 
223  // FIELD_B *= -0.5
224  // CHKERR fb->fieldScale(-0.5, "FIELD_B");
225  CHKERR fb->fieldLambdaOnValues(field_scale, "FIELD_B", &meshset_ents);
226 
227  // FIELD_B(i) *= FIELD_C
228  CHKERR fb->fieldLambdaOnEntities(vector_times_scalar_field, "FIELD_C",
229  "FIELD_B", true, &meshset_ents);
230 
231  auto dofs_ptr = m_field.get_dofs();
232  constexpr double eps = std::numeric_limits<double>::epsilon();
233  for (auto dit : *dofs_ptr) {
234  auto check = [&](const std::string name, const double expected) {
236  if (dit->getName() == name && dit->getEntType() == MBVERTEX) {
237  cout << name << " " << dit->getFieldData() << " " << expected
238  << endl;
239  if (abs(dit->getFieldData() - expected) > eps)
240  SETERRQ3(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
241  "Wrong DOF value %4.3e != %4.3e for %s",
242  dit->getFieldData(), expected,
243  boost::lexical_cast<std::string>(*dit).c_str());
244  }
245 
247  };
248 
249  CHKERR check("FIELD_A", 0.5);
250 
251  if (meshset_ents.find(dit->getEnt()) != meshset_ents.end()) {
252  CHKERR check("FIELD_B", 0.75);
253  CHKERR check("FIELD_C", 1.5);
254  }
255  }
256 
258  };
259 
260  CHKERR check_axpy();
261  CHKERR check_set_vertex();
262  CHKERR check_lambda();
263  }
264  CATCH_ERRORS;
265 
267  return 0;
268 }
help
static char help[]
Definition: field_blas_axpy.cpp:14
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
SIDESET
@ SIDESET
Definition: definitions.h:147
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Types::VectorDouble3
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::EntityMethod
Data structure to exchange data between mofem and User Loop Methods on entities.
Definition: LoopMethods.hpp:471
PRESSURESET
@ PRESSURESET
Definition: definitions.h:152
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Definition: MeshsetsManager.hpp:49
MoFEM.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::CoreInterface::get_dofs
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreTmp
Definition: Core.hpp:36
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
FTensor::dd
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
main
int main(int argc, char *argv[])
Definition: field_blas_axpy.cpp:16
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::CoreInterface::loop_entities
virtual MoFEMErrorCode loop_entities(const std::string field_name, EntityMethod &method, Range const *const ents=nullptr, int verb=DEFAULT_VERBOSITY)=0
Loop over field entities.