v0.14.0
body_force_atom_test.cpp
Go to the documentation of this file.
1 /** \file body_force_atom_test.cpp
2 
3  \brief Atom test for (linear) body forces implementation
4 
5 */
6 
7 
8 
10 using namespace MoFEM;
11 
12 namespace bio = boost::iostreams;
13 using bio::stream;
14 using bio::tee_device;
15 
16 static char help[] = "...\n\n";
17 
18 int main(int argc, char *argv[]) {
19 
20  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
21 
22  try {
23 
24  moab::Core mb_instance;
25  moab::Interface &moab = mb_instance;
26  int rank;
27  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
28 
29  PetscBool flg = PETSC_TRUE;
30  char mesh_file_name[255];
31  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
32  mesh_file_name, 255, &flg);
33  if (flg != PETSC_TRUE) {
34  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
35  }
36 
37  const char *option;
38  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
39  CHKERR moab.load_file(mesh_file_name, 0, option);
40 
41  // Create MoFEM (Joseph) database
42  MoFEM::Core core(moab);
43  MoFEM::Interface &m_field = core;
44 
45  // set entitities bit level
46  BitRefLevel bit_level0;
47  bit_level0.set(0);
48  EntityHandle meshset_level0;
49  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
50  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
51  0, 3, bit_level0);
52 
53  // Fields
54  CHKERR m_field.add_field("DISPLACEMENT", H1, AINSWORTH_LEGENDRE_BASE, 3);
55  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
56  3);
57 
58  // FE
59  CHKERR m_field.add_finite_element("TEST_FE");
60 
61  // Define rows/cols and element data
63  "DISPLACEMENT");
65  "DISPLACEMENT");
67  "DISPLACEMENT");
69  "MESH_NODE_POSITIONS");
70 
71  // Problem
72  CHKERR m_field.add_problem("TEST_PROBLEM");
73 
74  // set finite elements for problem
75  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
76  // set refinement level for problem
77  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
78 
79  // meshset consisting all entities in mesh
80  EntityHandle root_set = moab.get_root_set();
81  // add entities to field
82  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "DISPLACEMENT");
83  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
84  "MESH_NODE_POSITIONS");
85  // add entities to finite element
86  CHKERR m_field.add_ents_to_finite_element_by_type(root_set, MBTET,
87  "TEST_FE");
88 
89  // set app. order
90  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
91  // (Mark Ainsworth & Joe Coyle)
92  int order = 2;
93  CHKERR m_field.set_field_order(root_set, MBTET, "DISPLACEMENT", order);
94  CHKERR m_field.set_field_order(root_set, MBTRI, "DISPLACEMENT", order);
95  CHKERR m_field.set_field_order(root_set, MBEDGE, "DISPLACEMENT", order);
96  CHKERR m_field.set_field_order(root_set, MBVERTEX, "DISPLACEMENT", 1);
97 
98  CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
99  CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
100  CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
101  CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
102  1);
103 
104  /****/
105  // build database
106  // build field
107  CHKERR m_field.build_fields();
108  // build finite elemnts
109  CHKERR m_field.build_finite_elements();
110  // build adjacencies
111  CHKERR m_field.build_adjacencies(bit_level0);
112 
113  // build problem
114  ProblemsManager *prb_mng_ptr;
115  CHKERR m_field.getInterface(prb_mng_ptr);
116  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
117 
118  /****/
119  // mesh partitioning
120  // partition
121  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
122  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
123  // what are ghost nodes, see Petsc Manual
124  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
125 
126  // set from positions of 10 node tets
127  Projection10NodeCoordsOnField ent_method(m_field, "MESH_NODE_POSITIONS");
128  CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method);
129 
130  Vec F;
131  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
132  ROW, &F);
133 
134  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
135  typedef stream<TeeDevice> TeeStream;
136  std::ofstream ofs("body_force_atom_test.txt");
137  TeeDevice my_tee(std::cout, ofs);
138  TeeStream my_split(my_tee);
139 
140  CHKERR VecZeroEntries(F);
141  BodyForceConstantField body_forces_methods(m_field);
142 
144  m_field, BLOCKSET | BODYFORCESSET, it)) {
145  Block_BodyForces mydata;
146  CHKERR it->getAttributeDataStructure(mydata);
147  my_split << mydata << std::endl;
148  CHKERR body_forces_methods.addBlock("DISPLACEMENT", F,
149  it->getMeshsetId());
150  }
151  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE",
152  body_forces_methods.getLoopFe());
153  CHKERR VecAssemblyBegin(F);
154  CHKERR VecAssemblyEnd(F);
155 
156  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
157  "TEST_PROBLEM", COL, F, INSERT_VALUES, SCATTER_REVERSE);
158 
159  const Problem *problemPtr;
160  CHKERR m_field.get_problem("TEST_PROBLEM", &problemPtr);
161  std::map<EntityHandle, double> m0, m1, m2;
162  for (_IT_NUMEREDDOF_ROW_FOR_LOOP_(problemPtr, dit)) {
163 
164  if (dit->get()->getDofCoeffIdx() != 1)
165  continue;
166 
167  my_split.precision(3);
168  my_split.setf(std::ios::fixed);
169  my_split << dit->get()->getPetscGlobalDofIdx() << " "
170  << dit->get()->getFieldData() << std::endl;
171  }
172 
173  double sum = 0;
174  CHKERR VecSum(F, &sum);
175  my_split << std::endl
176  << "Sum : " << std::setprecision(3) << sum << std::endl;
177 
178  CHKERR VecDestroy(&F);
179  }
180  CATCH_ERRORS;
181 
183 
184  return 0;
185 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::Block_BodyForces
Body force data structure.
Definition: MaterialBlocks.hpp:313
TeeStream
stream< TeeDevice > TeeStream
Definition: forces_and_sources_testing_contact_prism_element.cpp:10
MoFEM::CoreInterface::loop_finite_elements
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.
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::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
BodyForceConstantField::addBlock
MoFEMErrorCode addBlock(const std::string field_name, Vec F, int ms_id)
Definition: BodyForce.hpp:94
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
_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::CoreInterface::get_problem
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
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
BodyForceConstantField::getLoopFe
MyVolumeFE & getLoopFe()
Definition: BodyForce.hpp:23
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
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::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
BodyForceConstantField
Body forces elements.
Definition: BodyForce.hpp:12
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
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.
BODYFORCESSET
@ BODYFORCESSET
block name is "BODY_FORCES"
Definition: definitions.h:162
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::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
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
_IT_NUMEREDDOF_ROW_FOR_LOOP_
#define _IT_NUMEREDDOF_ROW_FOR_LOOP_(PROBLEMPTR, IT)
use with loops to iterate row DOFs
Definition: ProblemsMultiIndices.hpp:164
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
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
BLOCKSET
@ BLOCKSET
Definition: definitions.h:148
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
help
static char help[]
Definition: body_force_atom_test.cpp:16
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
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::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
main
int main(int argc, char *argv[])
Definition: body_force_atom_test.cpp:18
F
@ F
Definition: free_surface.cpp:394