v0.14.0
Loading...
Searching...
No Matches
Functions | Variables
forces_and_sources_testing_vertex_element.cpp File Reference
#include <MoFEM.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 13 of file forces_and_sources_testing_vertex_element.cpp.

13 {
14
15 MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
16
17 try {
18
19 moab::Core mb_instance;
20 moab::Interface &moab = mb_instance;
21 int rank;
22 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23
24 PetscBool flg = PETSC_TRUE;
25 char mesh_file_name[255];
26#if PETSC_VERSION_GE(3, 6, 4)
27 CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
28 255, &flg);
29#else
30 CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
31 mesh_file_name, 255, &flg);
32#endif
33 if (flg != PETSC_TRUE) {
34 SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
35 }
36
37 // Create MoFEM (Joseph) database
38 MoFEM::Core core(moab);
39 MoFEM::Interface &m_field = core;
40
41 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
42 auto moab_comm_wrap =
43 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD, false);
44 if (pcomm == NULL)
45 pcomm =
46 new ParallelComm(&moab, moab_comm_wrap->get_comm());
47
48 const char *option;
49 option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
51 CHKERR moab.load_file(mesh_file_name, 0, option);
53
54 // set entitities bit level
55 BitRefLevel bit_level0;
56 bit_level0.set(0);
57 EntityHandle meshset_level0;
58 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
59 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
60 0, 3, bit_level0);
61
62 // Fields
63 CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 3);
64
65 // FE
66 CHKERR m_field.add_finite_element("TEST_FE");
67
68 // Define rows/cols and element data
69 CHKERR m_field.modify_finite_element_add_field_row("TEST_FE", "FIELD1");
70 CHKERR m_field.modify_finite_element_add_field_col("TEST_FE", "FIELD1");
71 CHKERR m_field.modify_finite_element_add_field_data("TEST_FE", "FIELD1");
72
73 // Problem
74 CHKERR m_field.add_problem("TEST_PROBLEM");
75
76 // set finite elements for problem
77 CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
78 // set refinement level for problem
79 CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
80
81 // meshset consisting all entities in mesh
82 EntityHandle root_set = moab.get_root_set();
83 // add entities to field
84 CHKERR m_field.add_ents_to_field_by_type(root_set, MBTET, "FIELD1");
85 // add entities to finite element
86 Range tets;
87 CHKERR moab.get_entities_by_type(0, MBTET, tets, false);
88 Skinner skin(&m_field.get_moab());
89 Range tets_skin;
90 CHKERR skin.find_skin(0, tets, false, tets_skin);
91 Range tets_skin_nodes;
92 CHKERR moab.get_connectivity(tets_skin, tets_skin_nodes, true);
93 CHKERR m_field.add_ents_to_finite_element_by_type(tets_skin_nodes, MBVERTEX,
94 "TEST_FE");
95
96 // set app. order
97 // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
98 // (Mark Ainsworth & Joe Coyle)
99 int order = 3;
100 CHKERR m_field.set_field_order(root_set, MBTET, "FIELD1", order);
101 CHKERR m_field.set_field_order(root_set, MBTRI, "FIELD1", order);
102 CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
103 CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
104
105 /****/
106 // build database
107 // build field
108 CHKERR m_field.build_fields();
109 // set FIELD1 from positions of 10 node tets
110 Projection10NodeCoordsOnField ent_method(m_field, "FIELD1");
111 CHKERR m_field.loop_dofs("FIELD1", ent_method);
112 // build finite elemnts
114 // build adjacencies
115 CHKERR m_field.build_adjacencies(bit_level0);
116 // build problem
117 ProblemsManager *prb_mng_ptr;
118 CHKERR m_field.getInterface(prb_mng_ptr);
119 CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
120
121 /****/
122 // mesh partitioning
123 // partition
124 CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
125 CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
126 // what are ghost nodes, see Petsc Manual
127 CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
128
130
131 typedef tee_device<std::ostream, std::ofstream> TeeDevice;
132 typedef stream<TeeDevice> TeeStream;
133
134 std::ofstream ofs("forces_and_sources_testing_vertex_element.txt");
135 TeeDevice my_tee(std::cout, ofs);
136 TeeStream my_split(my_tee);
137
139
140 TeeStream &my_split;
141 MyOp(TeeStream &_my_split, const char type)
143 "FIELD1", type),
144 my_split(_my_split) {}
145
146 MoFEMErrorCode doWork(int side, EntityType type,
149
150 my_split << "NH1" << std::endl;
151 my_split << "side: " << side << " type: " << type << std::endl;
152 my_split << "data: " << data << std::endl;
153 my_split << "coords: " << std::setprecision(3) << getCoords()
154 << std::endl;
155
157 }
158
159 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
160 EntityType col_type,
162 EntitiesFieldData::EntData &col_data) {
164 my_split << "ROW NH1NH1" << std::endl;
165 my_split << "row side: " << row_side << " row_type: " << row_type
166 << std::endl;
167 my_split << row_data << std::endl;
168 my_split << "COL NH1NH1" << std::endl;
169 my_split << "col side: " << col_side << " col_type: " << col_type
170 << std::endl;
171 my_split << col_data << std::endl;
173 }
174 };
175
176 fe1.getOpPtrVector().push_back(
177 new MyOp(my_split, ForcesAndSourcesCore::UserDataOperator::OPROW));
178 fe1.getOpPtrVector().push_back(
179 new MyOp(my_split, ForcesAndSourcesCore::UserDataOperator::OPROWCOL));
180
181 CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
182 }
184
186
187 return 0;
188}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define BARRIER_PCOMM_RANK_END(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MYPCOMM_INDEX
default communicator number PCOMM
#define BARRIER_PCOMM_RANK_START(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
tee_device< std::ostream, std::ofstream > TeeDevice
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 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 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_row(const std::string &fe_name, const std::string name_row)=0
set field row 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_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.
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 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
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
char mesh_file_name[255]
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
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:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:112
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MyOp(std::array< double, 12 > &eval_points)

Variable Documentation

◆ help

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