v0.14.0
field_approximation.cpp
Go to the documentation of this file.
1 
2 
4 
5 namespace bio = boost::iostreams;
6 using bio::stream;
7 using bio::tee_device;
8 
9 #define HOON
10 
11 static char help[] = "...\n\n";
12 
13 /// Example approx. function
14 struct MyFunApprox {
15 
16  std::vector<VectorDouble> result;
17 
18  std::vector<VectorDouble> &operator()(double x, double y, double z) {
19  result.resize(1);
20  result[0].resize(3);
21  (result[0])[0] = x;
22  (result[0])[1] = y;
23  (result[0])[2] = z * z;
24  return result;
25  }
26 };
27 
28 int main(int argc, char *argv[]) {
29 
30  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
31 
32  try {
33 
34  moab::Core mb_instance;
35  moab::Interface &moab = mb_instance;
36  int rank;
37  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
38 
39  PetscBool flg = PETSC_TRUE;
40  char mesh_file_name[255];
41  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
42  mesh_file_name, 255, &flg);
43  if (flg != PETSC_TRUE) {
44  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
45  "*** ERROR -my_file (MESH FILE NEEDED)");
46  }
47 
48  const char *option;
49  option = "";
50  CHKERR moab.load_file(mesh_file_name, 0, option);
51 
52  // Create MoFEM (Joseph) database
53  MoFEM::Core core(moab);
54  MoFEM::Interface &m_field = core;
55 
56  // set entities bit level
57  BitRefLevel bit_level0;
58  bit_level0.set(0);
59  EntityHandle meshset_level0;
60  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
61  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
62  0, 3, bit_level0);
63 
64  // Fields
65  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
66 #ifdef HOON
67  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
68  3, MB_TAG_SPARSE, MF_ZERO);
69 #endif
70 
71  // FE
72  CHKERR m_field.add_finite_element("TEST_FE");
73 
74  // Define rows/cols and element data
75  CHKERR m_field.modify_finite_element_add_field_row("TEST_FE", "FIELD1");
76  CHKERR m_field.modify_finite_element_add_field_col("TEST_FE", "FIELD1");
77  CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD1");
78 #ifdef HOON
80  "MESH_NODE_POSITIONS");
81 #endif
82 
83  // Problem
84  CHKERR m_field.add_problem("TEST_PROBLEM");
85 
86  // set finite elements for problem
87  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
88  // set refinement level for problem
89  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
90 
91  // meshset consisting all entities in mesh
92  EntityHandle root_set = moab.get_root_set();
93  // add entities to field
94  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
95 #ifdef HOON
96  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
97  "MESH_NODE_POSITIONS");
98 #endif
99  // add entities to finite element
100  CHKERR
101  m_field.add_ents_to_finite_element_by_type(root_set, MBTET, "TEST_FE");
102 
103  // set app. order
104  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
105  // (Mark Ainsworth & Joe Coyle)
106  int order = 3;
107  CHKERR
108  PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order, &flg);
109  if (flg != PETSC_TRUE) {
110  order = 3;
111  }
112  CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
113  CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
114  CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
115  CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
116 #ifdef HOON
117  CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
118  CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
119  CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
120  CHKERR
121  m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS", 1);
122 #endif
123 
124  /****/
125  // build database
126  // build field
127  CHKERR m_field.build_fields();
128 #ifdef HOON
129  Projection10NodeCoordsOnField ent_method_material(m_field,
130  "MESH_NODE_POSITIONS");
131  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method_material);
132 #endif
133  // build finite elemnts
134  CHKERR m_field.build_finite_elements();
135  // build adjacencies
136  CHKERR m_field.build_adjacencies(bit_level0);
137  // build problem
138  ProblemsManager *prb_mng_ptr;
139  CHKERR m_field.getInterface(prb_mng_ptr);
140  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
141 
142  /****/
143  // mesh partitioning
144  // partition
145  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
146  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
147  // what are ghost nodes, see Petsc Manual
148  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
149 
150  Mat A;
151  CHKERR
152  m_field.getInterface<MatrixManager>()
153  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("TEST_PROBLEM", &A);
154  Vec D, F;
155  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
156  ROW, &F);
157  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
158  COL, &D);
159 
160  std::vector<Vec> vec_F;
161  vec_F.push_back(F);
162 
163  {
164  MyFunApprox function_evaluator;
165  FieldApproximationH1 field_approximation(m_field);
166  field_approximation.loopMatrixAndVectorVolume(
167  "TEST_PROBLEM", "TEST_FE", "FIELD1", A, vec_F, function_evaluator);
168  }
169 
170  CHKERR MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
171  CHKERR MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
172  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
173  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
174  CHKERR VecGhostUpdateBegin(F, INSERT_VALUES, SCATTER_FORWARD);
175  CHKERR VecGhostUpdateEnd(F, INSERT_VALUES, SCATTER_FORWARD);
176 
177  KSP solver;
178  CHKERR KSPCreate(PETSC_COMM_WORLD, &solver);
179  CHKERR KSPSetOperators(solver, A, A);
180 
181  CHKERR PetscOptionsInsertString(NULL,
182  "-ksp_monitor -ksp_type fgmres -pc_type "
183  "bjacobi -ksp_atol 0 -ksp_rtol 1e-12");
184 
185  CHKERR KSPSetFromOptions(solver);
186  CHKERR PetscOptionsView(NULL, PETSC_VIEWER_STDOUT_WORLD);
187  CHKERR KSPSetUp(solver);
188  CHKERR KSPSolve(solver, F, D);
189  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
190  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
191 
192  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
193  "TEST_PROBLEM", COL, D, INSERT_VALUES, SCATTER_REVERSE);
194 
195  CHKERR KSPDestroy(&solver);
196  CHKERR VecDestroy(&D);
197  CHKERR VecDestroy(&F);
198  CHKERR MatDestroy(&A);
199 
200  EntityHandle fe_meshset = m_field.get_finite_element_meshset("TEST_FE");
201  Range tets;
202  CHKERR moab.get_entities_by_type(fe_meshset, MBTET, tets, true);
203  Range tets_edges;
204  CHKERR moab.get_adjacencies(tets, 1, false, tets_edges,
205  moab::Interface::UNION);
206  EntityHandle edges_meshset;
207  CHKERR moab.create_meshset(MESHSET_SET, edges_meshset);
208  CHKERR moab.add_entities(edges_meshset, tets);
209  CHKERR moab.add_entities(edges_meshset, tets_edges);
210  CHKERR moab.convert_entities(edges_meshset, true, false, false);
211 
212  ProjectionFieldOn10NodeTet ent_method_field1_on_10nodeTet(
213  m_field, "FIELD1", true, false, "FIELD1");
214  CHKERR m_field.loop_dofs("FIELD1", ent_method_field1_on_10nodeTet);
215  ent_method_field1_on_10nodeTet.setNodes = false;
216  CHKERR m_field.loop_dofs("FIELD1", ent_method_field1_on_10nodeTet);
217 
218  // if (pcomm->rank() == 0) {
219  // EntityHandle out_meshset;
220  // CHKERR moab.create_meshset(MESHSET_SET, out_meshset);
221  // CHKERR m_field.get_problem_finite_elements_entities(
222  // "TEST_PROBLEM", "TEST_FE", out_meshset);
223  // CHKERR moab.write_file("out.vtk", "VTK", "", &out_meshset, 1);
224  // CHKERR moab.delete_entities(&out_meshset, 1);
225  // }
226 
227  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
228  typedef stream<TeeDevice> TeeStream;
229 
230  std::ofstream ofs("field_approximation.txt");
231  TeeDevice tee(cout, ofs);
232  TeeStream my_split(tee);
233 
234  Range nodes;
235  CHKERR moab.get_entities_by_type(0, MBVERTEX, nodes, true);
236  MatrixDouble nodes_vals;
237  nodes_vals.resize(nodes.size(), 3);
238  CHKERR moab.tag_get_data(ent_method_field1_on_10nodeTet.th, nodes,
239  &*nodes_vals.data().begin());
240 
241  const double eps = 1e-4;
242 
243  my_split.precision(3);
244  my_split.setf(std::ios::fixed);
245  for (DoubleAllocator::iterator it = nodes_vals.data().begin();
246  it != nodes_vals.data().end(); it++) {
247  *it = fabs(*it) < eps ? 0.0 : *it;
248  }
249  my_split << nodes_vals << std::endl;
250 
251  const Problem *problemPtr;
252  CHKERR m_field.get_problem("TEST_PROBLEM", &problemPtr);
253  std::map<EntityHandle, double> m0, m1, m2;
254  for (_IT_NUMEREDDOF_ROW_FOR_LOOP_(problemPtr, dit)) {
255 
256  my_split.precision(3);
257  my_split.setf(std::ios::fixed);
258  double val = fabs(dit->get()->getFieldData()) < eps
259  ? 0.0
260  : dit->get()->getFieldData();
261  my_split << dit->get()->getPetscGlobalDofIdx() << " " << val << std::endl;
262  }
263  }
264  CATCH_ERRORS;
265 
267 
268  return 0;
269 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
TeeStream
stream< TeeDevice > TeeStream
Definition: forces_and_sources_testing_contact_prism_element.cpp:10
MoFEM::CoreInterface::loop_dofs
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreInterface::get_finite_element_meshset
virtual EntityHandle get_finite_element_meshset(const std::string name) const =0
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
FieldApproximationH1
Finite element for approximating analytical filed on the mesh.
Definition: FieldApproximation.hpp:17
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
BasicFiniteElements.hpp
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
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
TeeDevice
tee_device< std::ostream, std::ofstream > TeeDevice
Definition: forces_and_sources_testing_contact_prism_element.cpp:9
ROW
@ ROW
Definition: definitions.h:123
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
MyFunApprox
Example approx. function.
Definition: field_approximation.cpp:14
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::ProjectionFieldOn10NodeTet
Definition: Projection10NodeCoordsOnField.hpp:49
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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.
COL
@ COL
Definition: definitions.h:123
MoFEM::ProjectionFieldOn10NodeTet::setNodes
bool setNodes
Definition: Projection10NodeCoordsOnField.hpp:61
help
static char help[]
Definition: field_approximation.cpp:11
MyFunApprox::result
std::vector< VectorDouble > result
Definition: field_approximation.cpp:16
MoFEM::ProblemsManager::partitionFiniteElements
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
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
FieldApproximationH1::loopMatrixAndVectorVolume
MoFEMErrorCode loopMatrixAndVectorVolume(const std::string &problem_name, const std::string &fe_name, const std::string &field_name, Mat A, std::vector< Vec > &vec_F, FUNEVAL &function_evaluator)
assemble matrix and vector
Definition: FieldApproximation.hpp:481
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
MyFunApprox::operator()
std::vector< VectorDouble > & operator()(double x, double y, double z)
Definition: field_approximation.cpp:18
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
_IT_NUMEREDDOF_ROW_FOR_LOOP_
#define _IT_NUMEREDDOF_ROW_FOR_LOOP_(PROBLEMPTR, IT)
use with loops to iterate row DOFs
Definition: ProblemsMultiIndices.hpp:164
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
MF_ZERO
@ MF_ZERO
Definition: definitions.h:98
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
main
int main(int argc, char *argv[])
Definition: field_approximation.cpp:28
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
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
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
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::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
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
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::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
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::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::ProjectionFieldOn10NodeTet::th
Tag th
Definition: Projection10NodeCoordsOnField.hpp:67
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
F
@ F
Definition: free_surface.cpp:394