v0.13.0
prism_polynomial_approximation.cpp
Go to the documentation of this file.
1 /** \file prism_polynomial_approximation.cpp
2  \example prism_polynomial_approximation.cpp
3  \brief Checking approximation functions on prism
4 */
5 
6 /* This file is part of MoFEM.
7  * MoFEM is free software: you can redistribute it and/or modify it under
8  * the terms of the GNU Lesser General Public License as published by the
9  * Free Software Foundation, either version 3 of the License, or (at your
10  * option) any later version.
11  *
12  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
13  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  * License for more details.
16  *
17  * You should have received a copy of the GNU Lesser General Public
18  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
19 
20 #include <MoFEM.hpp>
21 
22 using namespace MoFEM;
23 
24 static char help[] = "...\n\n";
25 
26 static constexpr int approx_order = 6;
27 
29  static inline double fun(double x, double y, double z) {
30  double r = 1;
31  for (int o = 1; o <= approx_order; ++o) {
32  for (int i = 0; i <= o; ++i) {
33  for (int j = 0; j <= (o - i); ++j) {
34  int k = o - i - j;
35  if (k >= 0) {
36  r += pow(x, i) * pow(y, j) * pow(z, k);
37  }
38  }
39  }
40  }
41  return r;
42  }
43 
44  static inline VectorDouble3 diff_fun(double x, double y, double z) {
45  VectorDouble3 r(3);
46  r.clear();
47  for (int o = 1; o <= approx_order; ++o) {
48  for (int i = 0; i <= o; ++i) {
49  for (int j = 0; j <= (o - i); ++j) {
50  int k = o - i - j;
51  if (k >= 0) {
52  r[0] += i > 0 ? i * pow(x, i - 1) * pow(y, j) * pow(z, k) : 0;
53  r[1] += j > 0 ? j * pow(x, i) * pow(y, j - 1) * pow(z, k) : 0;
54  r[2] += k > 0 ? k * pow(x, i) * pow(y, j) * pow(z, k - 1) : 0;
55  }
56  }
57  }
58  }
59  return r;
60  }
61 };
62 
65 
66  PrismOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
67  boost::shared_ptr<MatrixDouble> &diff_field_vals);
68  MoFEMErrorCode doWork(int side, EntityType type,
70 
71 private:
72  boost::shared_ptr<VectorDouble> fieldVals;
73  boost::shared_ptr<MatrixDouble> diffFieldVals;
74 };
75 
76 struct PrismOpRhs
78 
80  MoFEMErrorCode doWork(int side, EntityType type,
82 
83 private:
85 };
86 
87 struct PrismOpLhs
89 
91  MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
92  EntityType col_type,
94  EntitiesFieldData::EntData &col_data);
95 
96 private:
98 };
99 
101 
106 };
107 
108 int main(int argc, char *argv[]) {
109 
110  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
111 
112  try {
113 
114  moab::Core mb_instance;
115  moab::Interface &moab = mb_instance;
116 
117  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
118  auto moab_comm_wrap =
119  boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
120  if (pcomm == NULL)
121  pcomm = new ParallelComm(&moab, moab_comm_wrap->get_comm());
122 
123  std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
124  0, 0, 1, 1, 0, 1, 0, 1, 1};
125  std::array<EntityHandle, 6> one_prism_nodes;
126  for (int n = 0; n != 6; ++n)
127  CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
128  EntityHandle one_prism;
129  CHKERR moab.create_element(MBPRISM, one_prism_nodes.data(), 6, one_prism);
130  Range one_prism_range;
131  one_prism_range.insert(one_prism);
132  Range one_prism_adj_ents;
133  for (int d = 1; d != 3; ++d)
134  CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
135  moab::Interface::UNION);
136 
137  MoFEM::Core core(moab);
138  MoFEM::Interface &m_field = core;
139 
140  PrismsFromSurfaceInterface *prisms_from_surface_interface;
141  CHKERR m_field.getInterface(prisms_from_surface_interface);
142  BitRefLevel bit_level0;
143  bit_level0.set(0);
144  CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
145  one_prism_range, bit_level0);
146  CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
147  bit_level0);
148 
149  // Fields
150  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
151  CHKERR m_field.add_ents_to_field_by_type(0, MBPRISM, "FIELD1");
152 
153  CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD1", 1);
154  CHKERR m_field.set_field_order(0, MBEDGE, "FIELD1", approx_order);
155  CHKERR m_field.set_field_order(0, MBTRI, "FIELD1", approx_order);
156  CHKERR m_field.set_field_order(0, MBQUAD, "FIELD1", approx_order);
157  CHKERR m_field.set_field_order(0, MBPRISM, "FIELD1", approx_order);
158  CHKERR m_field.build_fields();
159 
160  // FE
161  CHKERR m_field.add_finite_element("PRISM");
162 
163  // Define rows/cols and element data
164  CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
165  CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
166  CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
167  CHKERR m_field.add_ents_to_finite_element_by_type(0, MBPRISM, "PRISM");
168 
169  // build finite elemnts
170  CHKERR m_field.build_finite_elements();
171  // //build adjacencies
172  CHKERR m_field.build_adjacencies(bit_level0);
173 
174  // Problem
175  CHKERR m_field.add_problem("TEST_PROBLEM");
176 
177  // set finite elements for problem
178  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "PRISM");
179  // set refinement level for problem
180  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
181 
182  // build problem
183  ProblemsManager *prb_mng_ptr;
184  CHKERR m_field.getInterface(prb_mng_ptr);
185  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
186  // partition
187  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
188  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
189  // what are ghost nodes, see Petsc Manual
190  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
191 
192  // Create matrices
195  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", A);
197  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
198  ROW, F);
200  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
201  COL, D);
202 
203  auto assemble_matrices_and_vectors = [&]() {
205  PrismFE fe(m_field);
206  MatrixDouble inv_jac;
207  fe.getOpPtrVector().push_back(
209  fe.getOpPtrVector().push_back(
211  fe.getOpPtrVector().push_back(
212  new MoFEM::OpSetInvJacH1ForFatPrism(inv_jac));
213  fe.getOpPtrVector().push_back(new PrismOpRhs(F));
214  fe.getOpPtrVector().push_back(new PrismOpLhs(A));
215  CHKERR VecZeroEntries(F);
216  CHKERR MatZeroEntries(A);
217  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe);
218  CHKERR VecAssemblyBegin(F);
219  CHKERR VecAssemblyEnd(F);
220  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
221  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
223  };
224 
225  auto solve_problem = [&] {
227  auto solver = createKSP(PETSC_COMM_WORLD);
228  CHKERR KSPSetOperators(solver, A, A);
229  CHKERR KSPSetFromOptions(solver);
230  CHKERR KSPSetUp(solver);
231  CHKERR KSPSolve(solver, F, D);
232  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
233  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
234  CHKERR m_field.getInterface<VecManager>()->setLocalGhostVector(
235  "TEST_PROBLEM", COL, D, INSERT_VALUES, SCATTER_REVERSE);
237  };
238 
239  auto check_solution = [&] {
241  PrismFE fe(m_field);
242  boost::shared_ptr<VectorDouble> field_vals_ptr(new VectorDouble());
243  boost::shared_ptr<MatrixDouble> diff_field_vals_ptr(new MatrixDouble());
244  MatrixDouble inv_jac;
245  fe.getOpPtrVector().push_back(
246  new OpCalculateInvJacForFatPrism(inv_jac));
247  fe.getOpPtrVector().push_back(
248  new OpSetInvJacH1ForFatPrism(inv_jac));
249  fe.getOpPtrVector().push_back(
250  new OpCalculateScalarFieldValues("FIELD1", field_vals_ptr));
251  fe.getOpPtrVector().push_back(
252  new OpCalculateScalarFieldGradient<3>("FIELD1", diff_field_vals_ptr));
253  fe.getOpPtrVector().push_back(
254  new PrismOpCheck(field_vals_ptr, diff_field_vals_ptr));
255  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe);
257  };
258 
259  CHKERR assemble_matrices_and_vectors();
260  CHKERR solve_problem();
261  CHKERR check_solution();
262  }
263  CATCH_ERRORS;
264 
266  return 0;
267 }
268 
269 PrismOpCheck::PrismOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
270  boost::shared_ptr<MatrixDouble> &diff_field_vals)
272  "FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
273  fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
274 
278  if (type == MBVERTEX) {
279  const int nb_gauss_pts = data.getN().size2();
280  auto t_coords = getFTensor1CoordsAtGaussPts();
281  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
282  double f = ApproxFunction::fun(t_coords(0), t_coords(1), t_coords(2));
283  VectorDouble3 diff_f =
284  ApproxFunction::diff_fun(t_coords(0), t_coords(1), t_coords(2));
285 
286  std::cout << f - (*fieldVals)[gg] << " : ";
287  for (auto d : {0, 1, 2})
288  std::cout << diff_f[d] - (*diffFieldVals)(d, gg) << " ";
289  std::cout << std::endl;
290 
291  constexpr double eps = 1e-6;
292  if (std::abs(f - (*fieldVals)[gg]) > eps ||
293  !std::isnormal((*fieldVals)[gg]))
294  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
295  "Wrong value %6.4e != %6.4e (%6.4e)", f, (*fieldVals)[gg],
296  f - (*fieldVals)[gg]);
297 
298  for (auto d : {0, 1, 2})
299  if (std::abs(diff_f[d] - (*diffFieldVals)(d, gg)) > eps ||
300  !std::isnormal((*diffFieldVals)(d, gg)))
301  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
302  "Wrong diff value %6.4e != %6.4e (%6.4e)", diff_f[d],
303  (*diffFieldVals)(d, gg),
304  diff_f[d] - (*diffFieldVals)(d, gg));
305 
306  ++t_coords;
307  }
308  }
310 }
311 
314  "FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
315  F(f) {}
316 
320  const int nb_dofs = data.getN().size2();
321  if (nb_dofs) {
322  const int nb_gauss_pts = data.getN().size1();
323  VectorDouble nf(nb_dofs);
324  nf.clear();
325  auto t_base = data.getFTensor0N();
326  auto t_coords = getFTensor1CoordsAtGaussPts();
327  auto t_w = getFTensor0IntegrationWeight();
328  double vol = getVolume();
329  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
330  double f = ApproxFunction::fun(t_coords(0), t_coords(1), t_coords(2));
331  double v = t_w * vol * f;
332  double *val = &*nf.begin();
333  for (int bb = 0; bb != nb_dofs; ++bb) {
334  *val += v * t_base;
335  ++t_base;
336  ++val;
337  }
338  ++t_coords;
339  ++t_w;
340  }
341  CHKERR VecSetValues(F, data, &*nf.data().begin(), ADD_VALUES);
342  }
344 }
345 
348  "FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROWCOL),
349  A(a) {
350  // FIXME: Can be symmetric, is not for simplicity
351  sYmm = false;
352 }
353 
354 MoFEMErrorCode PrismOpLhs::doWork(int row_side, int col_side,
355  EntityType row_type, EntityType col_type,
356  EntitiesFieldData::EntData &row_data,
357  EntitiesFieldData::EntData &col_data) {
359  const int row_nb_dofs = row_data.getN().size2();
360  const int col_nb_dofs = col_data.getN().size2();
361  if (row_nb_dofs && col_nb_dofs) {
362  const int nb_gauss_pts = row_data.getN().size1();
363  MatrixDouble m(row_nb_dofs, col_nb_dofs);
364  m.clear();
365  auto t_w = getFTensor0IntegrationWeight();
366  double vol = getVolume();
367  double *row_base_ptr = &*row_data.getN().data().begin();
368  double *col_base_ptr = &*col_data.getN().data().begin();
369  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
370 
371  double v = t_w * vol;
372  cblas_dger(CblasRowMajor, row_nb_dofs, col_nb_dofs, v, row_base_ptr, 1,
373  col_base_ptr, 1, &*m.data().begin(), col_nb_dofs);
374 
375  row_base_ptr += row_nb_dofs;
376  col_base_ptr += col_nb_dofs;
377  ++t_w;
378  }
379  CHKERR MatSetValues(A, row_data, col_data, &*m.data().begin(), ADD_VALUES);
380  }
382 }
383 
384 int PrismFE::getRuleTrianglesOnly(int order) { return 2 * (order + 1); };
385 int PrismFE::getRuleThroughThickness(int order) { return 2 * (order + 1); };
static Index< 'n', 3 > n
static Index< 'o', 3 > o
ForcesAndSourcesCore::UserDataOperator UserDataOperator
constexpr double a
static const double eps
@ COL
Definition: definitions.h:136
@ ROW
Definition: definitions.h:136
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:73
@ H1
continuous field
Definition: definitions.h:98
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:53
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
virtual MoFEMErrorCode loop_finite_elements(const std::string &problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string &name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
FTensor::Index< 'i', SPACE_DIM > i
double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:103
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:88
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
UBlasVector< double > VectorDouble
Definition: Types.hpp:79
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
auto createKSP
CoreTmp< 0 > Core
Definition: Core.hpp:1096
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
const double D
diffusivity
const double r
rate factor
double A
int main(int argc, char *argv[])
static char help[]
static constexpr int approx_order
FTensor::Index< 'm', 3 > m
static VectorDouble3 diff_fun(double x, double y, double z)
static double fun(double x, double y, double z)
Managing BitRefLevels.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
Core (interface) class.
Definition: Core.hpp:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
structure to get information form mofem into EntitiesFieldData
boost::ptr_vector< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Matrix manager is used to build and partition problems.
Calculate inverse of jacobian for face element.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Operator for fat prism element updating integration weights in the volume.
Transform local reference derivatives of shape functions to global derivatives.
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:33
int getRuleThroughThickness(int order)
int getRuleTrianglesOnly(int order)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > diffFieldVals
boost::shared_ptr< VectorDouble > fieldVals
PrismOpCheck(boost::shared_ptr< VectorDouble > &field_vals, boost::shared_ptr< MatrixDouble > &diff_field_vals)
SmartPetscObj< Mat > A
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
PrismOpLhs(SmartPetscObj< Mat > &a)
SmartPetscObj< Vec > F
PrismOpRhs(SmartPetscObj< Vec > &f)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.