|
| v0.14.0
|
Go to the documentation of this file.
6 namespace bio = boost::iostreams;
10 static char help[] =
"...\n\n";
12 int main(
int argc,
char *argv[]) {
21 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
23 PetscBool flg = PETSC_TRUE;
27 if (flg != PETSC_TRUE) {
28 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
43 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
62 "MESH_NODE_POSITIONS");
67 std::ostringstream fe_name;
68 fe_name <<
"FORCE_FE_" << it->getMeshsetId();
80 CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris,
true);
88 std::ostringstream fe_name;
89 fe_name <<
"PRESSURE_FE_" << it->getMeshsetId();
98 fe_name.str(),
"MESH_NODE_POSITIONS");
103 CHKERR moab.get_entities_by_type(it->meshset, MBTRI, tris,
true);
151 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
153 std::ostringstream txt_name;
155 std::ofstream ofs(txt_name.str().c_str());
160 boost::ptr_map<std::string, NeumannForcesSurface> neumann_forces;
163 std::ostringstream fe_name;
164 fe_name <<
"FORCE_FE_" << it->getMeshsetId();
165 string fe_name_str = fe_name.str();
168 "MESH_NODE_POSITIONS", neumann_forces.at(fe_name_str).getLoopFe(),
170 neumann_forces.at(fe_name_str)
171 .addForce(
"DISPLACEMENT",
F, it->getMeshsetId());
173 CHKERR it->getBcDataStructure(data);
174 my_split << *it << std::endl;
175 my_split << data << std::endl;
179 std::ostringstream fe_name;
180 fe_name <<
"PRESSURE_FE_" << it->getMeshsetId();
181 string fe_name_str = fe_name.str();
184 "MESH_NODE_POSITIONS", neumann_forces.at(fe_name_str).getLoopFe(),
186 neumann_forces.at(fe_name_str)
187 .addPressure(
"DISPLACEMENT",
F, it->getMeshsetId());
189 CHKERR it->getBcDataStructure(data);
190 my_split << *it << std::endl;
191 my_split << data << std::endl;
193 boost::ptr_map<std::string, NeumannForcesSurface>::iterator mit =
194 neumann_forces.begin();
195 for (; mit != neumann_forces.end(); mit++) {
197 mit->second->getLoopFe());
203 "TEST_PROBLEM",
ROW,
F, INSERT_VALUES, SCATTER_REVERSE);
205 const double eps = 1e-4;
211 my_split.precision(3);
212 my_split.setf(std::ios::fixed);
213 double val = fabs(dit->get()->getFieldData()) <
eps
215 : dit->get()->getFieldData();
216 my_split << dit->get()->getPetscGlobalDofIdx() <<
" " << val << std::endl;
221 sum = fabs(sum) <
eps ? 0.0 : sum;
222 my_split << std::endl
223 <<
"Sum : " << std::setprecision(3) << sum << std::endl;
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
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.
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.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Problem manager is used to build and partition problems.
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
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
virtual const Problem * get_problem(const std::string problem_name) const =0
Get the problem object.
Definition of the force bc data structure.
Projection of edge entities with one mid-node on hierarchical basis.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
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.
Definition of the pressure bc data structure.
Deprecated interface functions.
DeprecatedCoreInterface Interface
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
#define CHKERR
Inline error check.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
implementation of Data Operators for Forces and Sources
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 build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
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 modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Vector manager is used to create vectors \mofem_vectors.
#define _IT_NUMEREDDOF_ROW_FOR_LOOP_(PROBLEMPTR, IT)
use with loops to iterate row DOFs
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
int main(int argc, char *argv[])
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
const FTensor::Tensor2< T, Dim, Dim > Vec
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
keeps basic data about problem
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.
Finite element and operators to apply force/pressures applied to surfaces.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
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.