v0.14.0
Functions | Variables
neumann_forces.cpp File Reference
#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 12 of file neumann_forces.cpp.

12  {
13 
14  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
15 
16  try {
17 
18  moab::Core mb_instance;
19  moab::Interface &moab = mb_instance;
20  int rank;
21  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
22 
23  PetscBool flg = PETSC_TRUE;
24  char mesh_file_name[255];
25  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
26  mesh_file_name, 255, &flg);
27  if (flg != PETSC_TRUE) {
28  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
29  }
30 
31  const char *option;
32  option = "";
33  CHKERR moab.load_file(mesh_file_name, 0, option);
34 
35  // Create MoFEM (Joseph) database
36  MoFEM::Core core(moab);
37  MoFEM::Interface &m_field = core;
38 
39  // set entitities bit level
40  BitRefLevel bit_level0;
41  bit_level0.set(0);
42  EntityHandle meshset_level0;
43  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
44  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
45  0, 3, bit_level0);
46 
47  // Fields
48  CHKERR m_field.add_field("DISPLACEMENT", H1, AINSWORTH_LEGENDRE_BASE, 3);
49  CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
50  3);
51 
52  // Problem
53  CHKERR m_field.add_problem("TEST_PROBLEM");
54  // set refinement level for problem
55  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
56 
57  // meshset consisting all entities in mesh
58  EntityHandle root_set = moab.get_root_set();
59  // add entities to field
60  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "DISPLACEMENT");
61  CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET,
62  "MESH_NODE_POSITIONS");
63 
65  it)) {
66 
67  std::ostringstream fe_name;
68  fe_name << "FORCE_FE_" << it->getMeshsetId();
69  CHKERR m_field.add_finite_element(fe_name.str());
70  CHKERR m_field.modify_finite_element_add_field_row(fe_name.str(),
71  "DISPLACEMENT");
72  CHKERR m_field.modify_finite_element_add_field_col(fe_name.str(),
73  "DISPLACEMENT");
74  CHKERR m_field.modify_finite_element_add_field_data(fe_name.str(),
75  "DISPLACEMENT");
76  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
77  fe_name.str());
78 
79  Range tris;
80  CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris, true);
81  CHKERR m_field.add_ents_to_finite_element_by_type(tris, MBTRI,
82  fe_name.str());
83  }
84 
86  m_field, SIDESET | PRESSURESET, it)) {
87 
88  std::ostringstream fe_name;
89  fe_name << "PRESSURE_FE_" << it->getMeshsetId();
90  CHKERR m_field.add_finite_element(fe_name.str());
91  CHKERR m_field.modify_finite_element_add_field_row(fe_name.str(),
92  "DISPLACEMENT");
93  CHKERR m_field.modify_finite_element_add_field_col(fe_name.str(),
94  "DISPLACEMENT");
95  CHKERR m_field.modify_finite_element_add_field_data(fe_name.str(),
96  "DISPLACEMENT");
98  fe_name.str(), "MESH_NODE_POSITIONS");
99  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM",
100  fe_name.str());
101 
102  Range tris;
103  CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris, true);
104  CHKERR m_field.add_ents_to_finite_element_by_type(tris, MBTRI,
105  fe_name.str());
106  }
107 
108  // set app. order
109  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
110  // (Mark Ainsworth & Joe Coyle)
111  int order = 2;
112  CHKERR m_field.set_field_order(root_set, MBTET, "DISPLACEMENT", order);
113  CHKERR m_field.set_field_order(root_set, MBTRI, "DISPLACEMENT", order);
114  CHKERR m_field.set_field_order(root_set, MBEDGE, "DISPLACEMENT", order);
115  CHKERR m_field.set_field_order(root_set, MBVERTEX, "DISPLACEMENT", 1);
116  CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 2);
117  CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 2);
118  CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 2);
119  CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
120  1);
121 
122  /****/
123  // build database
124  // build field
125  CHKERR m_field.build_fields();
126  // set FIELD1 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  // build finite elemnts
131  CHKERR m_field.build_finite_elements();
132  // build adjacencies
133  CHKERR m_field.build_adjacencies(bit_level0);
134  // build problem
135  ProblemsManager *prb_mng_ptr;
136  CHKERR m_field.getInterface(prb_mng_ptr);
137  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
138 
139  /****/
140  // mesh partitioning
141  // partition
142  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
143  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
144  // what are ghost nodes, see Petsc Manual
145  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
146 
147  Vec F;
148  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost("TEST_PROBLEM",
149  ROW, &F);
150 
151  typedef tee_device<std::ostream, std::ofstream> TeeDevice;
152  typedef stream<TeeDevice> TeeStream;
153  std::ostringstream txt_name;
154  txt_name << "forces_and_sources_" << mesh_file_name << ".txt";
155  std::ofstream ofs(txt_name.str().c_str());
156  TeeDevice my_tee(std::cout, ofs);
157  TeeStream my_split(my_tee);
158 
159  CHKERR VecZeroEntries(F);
160  boost::ptr_map<std::string, NeumannForcesSurface> neumann_forces;
162  it)) {
163  std::ostringstream fe_name;
164  fe_name << "FORCE_FE_" << it->getMeshsetId();
165  string fe_name_str = fe_name.str();
166  neumann_forces.insert(fe_name_str, new NeumannForcesSurface(m_field));
168  "MESH_NODE_POSITIONS", neumann_forces.at(fe_name_str).getLoopFe(),
169  false, false);
170  neumann_forces.at(fe_name_str)
171  .addForce("DISPLACEMENT", F, it->getMeshsetId());
172  ForceCubitBcData data;
173  CHKERR it->getBcDataStructure(data);
174  my_split << *it << std::endl;
175  my_split << data << std::endl;
176  }
178  m_field, SIDESET | PRESSURESET, it)) {
179  std::ostringstream fe_name;
180  fe_name << "PRESSURE_FE_" << it->getMeshsetId();
181  string fe_name_str = fe_name.str();
182  neumann_forces.insert(fe_name_str, new NeumannForcesSurface(m_field));
184  "MESH_NODE_POSITIONS", neumann_forces.at(fe_name_str).getLoopFe(),
185  false, false);
186  neumann_forces.at(fe_name_str)
187  .addPressure("DISPLACEMENT", F, it->getMeshsetId());
188  PressureCubitBcData data;
189  CHKERR it->getBcDataStructure(data);
190  my_split << *it << std::endl;
191  my_split << data << std::endl;
192  }
193  boost::ptr_map<std::string, NeumannForcesSurface>::iterator mit =
194  neumann_forces.begin();
195  for (; mit != neumann_forces.end(); mit++) {
196  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", mit->first,
197  mit->second->getLoopFe());
198  }
199  CHKERR VecAssemblyBegin(F);
200  CHKERR VecAssemblyEnd(F);
201 
202  CHKERR m_field.getInterface<VecManager>()->setGlobalGhostVector(
203  "TEST_PROBLEM", ROW, F, INSERT_VALUES, SCATTER_REVERSE);
204 
205  const double eps = 1e-4;
206 
207  const Problem *problemPtr;
208  CHKERR m_field.get_problem("TEST_PROBLEM", &problemPtr);
209  for (_IT_NUMEREDDOF_ROW_FOR_LOOP_(problemPtr, dit)) {
210 
211  my_split.precision(3);
212  my_split.setf(std::ios::fixed);
213  double val = fabs(dit->get()->getFieldData()) < eps
214  ? 0.0
215  : dit->get()->getFieldData();
216  my_split << dit->get()->getPetscGlobalDofIdx() << " " << val << std::endl;
217  }
218 
219  double sum = 0;
220  CHKERR VecSum(F, &sum);
221  sum = fabs(sum) < eps ? 0.0 : sum;
222  my_split << std::endl
223  << "Sum : " << std::setprecision(3) << sum << std::endl;
224 
225  CHKERR VecDestroy(&F);
226  }
227  CATCH_ERRORS;
228 
230 
231  return 0;
232 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 10 of file neumann_forces.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::addHOOpsFace3D
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
Definition: HODataOperators.hpp:699
SIDESET
@ SIDESET
Definition: definitions.h:160
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
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
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.
MoFEM::ForceCubitBcData
Definition of the force bc data structure.
Definition: BCData.hpp:139
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.
MoFEM::PressureCubitBcData
Definition of the pressure bc data structure.
Definition: BCData.hpp:379
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:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
NODESET
@ NODESET
Definition: definitions.h:159
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:548
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
FORCESET
@ FORCESET
Definition: definitions.h:164
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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_field)=0
set finite element field data
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:22
_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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
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::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.
NeumannForcesSurface
Finite element and operators to apply force/pressures applied to surfaces.
Definition: SurfacePressure.hpp:14
help
static char help[]
Definition: neumann_forces.cpp:10
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.
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.
F
@ F
Definition: free_surface.cpp:394