v0.14.0
Loading...
Searching...
No Matches
forces_and_sources_testing_edge_element.cpp

Testing edge elements nodesets.

Testing edge elements nodesets.

/** \file forces_and_sources_testing_edge_element.cpp
* \example forces_and_sources_testing_edge_element.cpp
* \brief Testing edge elements
* nodesets.
*
*/
#include <MoFEM.hpp>
namespace bio = boost::iostreams;
using bio::stream;
using bio::tee_device;
using namespace MoFEM;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
// Create MoFEM database
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
const char *option;
option = ""; ;
CHKERR moab.load_file(mesh_file_name, 0, option);
// set ebturues bit level
BitRefLevel bit_level0;
bit_level0.set(0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
enum bases { AINSWORTH, BERNSTEIN_BEZIER, LASBASETOP };
const char *list_bases[] = {"ainsworth", "bernstein_bezier"};
PetscInt choice_base_value = AINSWORTH;
CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
LASBASETOP, &choice_base_value, &flg);
if (choice_base_value == AINSWORTH)
else if (choice_base_value == BERNSTEIN_BEZIER)
// Fields
CHKERR m_field.add_field("FIELD1", H1, base, 3);
CHKERR m_field.add_field("MESH_NODE_POSITIONS", H1, AINSWORTH_LEGENDRE_BASE,
3);
// FE
CHKERR m_field.add_finite_element("TEST_FE");
// Define rows/cols and element data
auto add_field_to_fe = [&m_field](const std::string field_name) {
};
CHKERR add_field_to_fe("FIELD1");
"MESH_NODE_POSITIONS");
// Problem
CHKERR m_field.add_problem("TEST_PROBLEM");
// set finite elements for problem
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TEST_FE");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
// meshset consisting all entities in mesh
EntityHandle root_set = moab.get_root_set();
// add entities to field
CHKERR m_field.add_ents_to_field_by_type(root_set, MBEDGE, "FIELD1");
CHKERR m_field.add_ents_to_field_by_type(root_set, MBEDGE,
"MESH_NODE_POSITIONS");
CHKERR m_field.set_field_order(root_set, MBVERTEX, "MESH_NODE_POSITIONS",
1);
CHKERR m_field.set_field_order(root_set, MBEDGE, "MESH_NODE_POSITIONS", 1);
CHKERR m_field.set_field_order(root_set, MBTRI, "MESH_NODE_POSITIONS", 1);
CHKERR m_field.set_field_order(root_set, MBTET, "MESH_NODE_POSITIONS", 1);
// add entities to finite element
Range tets;
CHKERR moab.get_entities_by_type(0, MBTET, tets, false);
Skinner skin(&m_field.get_moab());
Range tets_skin;
CHKERR skin.find_skin(0, tets, false, tets_skin);
Range tets_skin_edges;
CHKERR moab.get_adjacencies(tets_skin, 1, false, tets_skin_edges,
moab::Interface::UNION);
CHKERR m_field.add_ents_to_finite_element_by_type(tets_skin_edges, MBEDGE,
"TEST_FE");
// set app. order
// see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
// (Mark Ainsworth & Joe Coyle)
int order = 3;
CHKERR m_field.set_field_order(root_set, MBEDGE, "FIELD1", order);
CHKERR m_field.set_field_order(root_set, MBVERTEX, "FIELD1", 1);
/****/
// build database
// build field
CHKERR m_field.build_fields();
// set FIELD1 from positions of 10 node tets
if (base == AINSWORTH_LEGENDRE_BASE) {
Projection10NodeCoordsOnField ent_method(m_field, "FIELD1");
CHKERR m_field.loop_dofs("FIELD1", ent_method);
}
{
Projection10NodeCoordsOnField ent_method(m_field, "MESH_NODE_POSITIONS");
CHKERR m_field.loop_dofs("MESH_NODE_POSITIONS", ent_method);
}
// build finite elemnts
// build adjacencies
CHKERR m_field.build_adjacencies(bit_level0);
// build problem
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
// mesh partitioning
// partition
CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
// what are ghost nodes, see Petsc Manual
CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
typedef tee_device<std::ostream, std::ofstream> TeeDevice;
typedef stream<TeeDevice> TeeStream;
std::ofstream ofs(
(std ::string("forces_and_sources_testing_edge_element_") +
ApproximationBaseNames[base] + ".txt")
.c_str());
TeeDevice my_tee(std::cout, ofs);
TeeStream my_split(my_tee);
struct MyOp : public EdgeElementForcesAndSourcesCore::UserDataOperator {
TeeStream &my_split;
MyOp(TeeStream &_my_split, const char type)
: EdgeElementForcesAndSourcesCore::UserDataOperator("FIELD1",
"FIELD1", type),
my_split(_my_split) {}
my_split << "NH1" << std::endl;
my_split << "side: " << side << " type: " << type << std::endl;
my_split << "data: " << data << std::endl;
my_split << "coords: " << std::setprecision(3) << getCoords()
<< std::endl;
my_split << "getCoordsAtGaussPts: " << std::setprecision(3)
<< getCoordsAtGaussPts() << std::endl;
my_split << "length: " << std::setprecision(3) << getLength()
<< std::endl;
my_split << "direction: " << std::setprecision(3) << getDirection()
<< std::endl;
int nb_gauss_pts = data.getN().size1();
for (int gg = 0; gg < nb_gauss_pts; gg++) {
my_split << "tangent " << gg << " " << getTangentAtGaussPts()
<< std::endl;
}
}
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
EntityType col_type,
my_split << "ROW NH1NH1" << std::endl;
my_split << "row side: " << row_side << " row_type: " << row_type
<< std::endl;
my_split << row_data << std::endl;
my_split << "COL NH1NH1" << std::endl;
my_split << "col side: " << col_side << " col_type: " << col_type
<< std::endl;
my_split << col_data << std::endl;
}
};
fe1.getOpPtrVector().push_back(
new OpGetHOTangentsOnEdge("MESH_NODE_POSITIONS"));
fe1.getOpPtrVector().push_back(
new MyOp(my_split, ForcesAndSourcesCore::UserDataOperator::OPROW));
fe1.getOpPtrVector().push_back(
new MyOp(my_split, ForcesAndSourcesCore::UserDataOperator::OPROWCOL));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TEST_FE", fe1);
}
return 0;
}
static char help[]
int main()
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
Definition definitions.h:58
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ NOBASE
Definition definitions.h:59
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition definitions.h:64
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
static const char *const ApproximationBaseNames[]
Definition definitions.h:72
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
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.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr auto field_name
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)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
Calculate tangent vector on edge form HO geometry approximation.
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)