static char help[] =
"teting interface inserting algorithm\n\n";
constexpr
bool debug =
false;
int main(
int argc,
char *argv[]) {
try {
PetscBool flg = PETSC_TRUE;
#if PETSC_VERSION_GE(3, 6, 4)
255, &flg);
#else
#endif
if (flg != PETSC_TRUE) {
"*** ERROR -my_file (MESH FILE NEEDED)");
}
const char *option;
option = "";
PrismInterface *interface;
std::vector<BitRefLevel> bit_levels;
int ll = 1;
CHKERR PetscPrintf(PETSC_COMM_WORLD,
"Insert Interface %d\n",
cit->getMeshsetId());
{
CHKERR moab.create_meshset(MESHSET_SET, ref_level_meshset);
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
ref_level_meshset);
->getEntitiesByTypeAndRefLevel(bit_levels.back(),
ref_level_meshset);
Range ref_level_tets;
CHKERR moab.get_entities_by_handle(ref_level_meshset, ref_level_tets,
true);
CHKERR interface->getSides(cubit_meshset, bit_levels.back(),
true, 0);
CHKERR interface->splitSides(ref_level_meshset, bit_levels.back(),
cubit_meshset, true, true, 0);
CHKERR moab.delete_entities(&ref_level_meshset, 1);
}
->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
cubit_meshset, MBVERTEX, true);
->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
cubit_meshset, MBEDGE, true);
->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
cubit_meshset, MBTRI, true);
->updateMeshsetByEntitiesChildren(cubit_meshset, bit_levels.back(),
cubit_meshset, MBTET, true);
}
}
"H1FIELD_SCALAR");
"H1FIELD_SCALAR");
"H1FIELD_SCALAR");
"H1FIELD_SCALAR");
"H1FIELD_SCALAR");
"H1FIELD_SCALAR");
true);
true);
for (int lll = ll - 2; lll < ll; lll++) {
std::stringstream problem_name;
problem_name << "PROBLEM_SCALAR_" << lll;
"ELEM_SCALAR");
"INTERFACE");
}
for (int lll = ll - 2; lll < ll; lll++) {
std::stringstream problem_name;
problem_name << "PROBLEM_SCALAR_" << lll;
std::stringstream message;
message << "set problem problem < " << problem_name.str()
<< " > bit level " << bit_levels[lll] << std::endl;
CHKERR PetscPrintf(PETSC_COMM_WORLD, message.str().c_str());
bit_levels[lll]);
}
Range tets_back_bit_level;
bit_levels.back(),
BitRefLevel().set(), tets_back_bit_level);
BlockSetAttributes mydata;
CHKERR cit->getAttributeDataStructure(mydata);
std::cout << mydata << std::endl;
Range tets;
CHKERR moab.get_entities_by_type(cubit_meshset, MBTET, tets,
true);
tets = intersect(tets_back_bit_level, tets);
Range nodes;
CHKERR moab.get_connectivity(tets, nodes,
true);
for (Range::iterator nit = nodes.begin(); nit != nodes.end(); nit++) {
double coords[3];
CHKERR moab.get_coords(&*nit, 1, coords);
coords[0] += mydata.data.User1;
coords[1] += mydata.data.User2;
coords[2] += mydata.data.User3;
CHKERR moab.set_coords(&*nit, 1, coords);
}
}
ProblemsManager *prb_mng_ptr;
for (int lll = ll - 2; lll < ll; lll++) {
std::stringstream problem_name;
problem_name << "PROBLEM_SCALAR_" << lll;
CHKERR prb_mng_ptr->buildProblem(problem_name.str(),
true);
CHKERR prb_mng_ptr->partitionProblem(problem_name.str());
CHKERR prb_mng_ptr->partitionFiniteElements(problem_name.str());
CHKERR prb_mng_ptr->partitionGhostDofs(problem_name.str());
}
std::ofstream myfile;
CHKERR moab.create_meshset(MESHSET_SET, out_meshset_tet);
bit_levels.back(),
BitRefLevel().set(), MBTET, out_meshset_tet);
Range tets;
CHKERR moab.get_entities_by_handle(out_meshset_tet, tets,
true);
for (Range::iterator tit = tets.begin(); tit != tets.end(); tit++) {
int num_nodes;
CHKERR moab.get_connectivity(*tit, conn, num_nodes,
true);
for (int nn = 0; nn < num_nodes; nn++) {
std::cout << conn[nn] << " ";
myfile << conn[nn] << " ";
}
std::cout << std::endl;
myfile << std::endl;
}
CHKERR moab.create_meshset(MESHSET_SET, out_meshset_prism);
bit_levels.back(),
BitRefLevel().set(), MBPRISM, out_meshset_prism);
Range prisms;
CHKERR moab.get_entities_by_handle(out_meshset_prism, prisms);
for (Range::iterator pit = prisms.begin(); pit != prisms.end(); pit++) {
int num_nodes;
CHKERR moab.get_connectivity(*pit, conn, num_nodes,
true);
for (int nn = 0; nn < num_nodes; nn++) {
std::cout << conn[nn] << " ";
myfile << conn[nn] << " ";
}
std::cout << std::endl;
myfile << std::endl;
}
myfile.close();
CHKERR moab.write_file(
"out_tet.vtk",
"VTK",
"", &out_meshset_tet, 1);
CHKERR moab.write_file(
"out_prism.vtk",
"VTK",
"", &out_meshset_prism, 1);
CHKERR moab.create_meshset(MESHSET_SET, out_meshset_tets_and_prism);
CHKERR moab.add_entities(out_meshset_tets_and_prism, tets);
CHKERR moab.add_entities(out_meshset_tets_and_prism, prisms);
CHKERR moab.write_file(
"out_tets_and_prisms.vtk",
"VTK",
"",
&out_meshset_tets_and_prism, 1);
CHKERR moab.create_meshset(MESHSET_SET, out_meshset_tris);
Range tris;
CHKERR moab.get_adjacencies(prisms, 2,
false, tris,
moab::Interface::UNION);
std::cerr << tris.size() << " : " << prisms.size() << std::endl;
CHKERR moab.add_entities(out_meshset_tris, tris);
CHKERR moab.write_file(
"out_tris.vtk",
"VTK",
"", &out_meshset_tris, 1);
}
}
}
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define CHKERR
Inline error check.
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 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 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_col(const std::string &fe_name, const std::string &name_row)=0
set field col 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.
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
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
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
int main(int argc, char *argv[])
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
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.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.