12static char help[] =
"Read file and print boundary conditions (ex. "
13 "./cubit_bc_test -my_file disp01.h5m) \n\n";
15int main(
int argc,
char *argv[]) {
21 moab::Core mb_instance;
22 moab::Interface &moab = mb_instance;
24 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
27 PetscBool flg = PETSC_TRUE;
29#if PETSC_VERSION_GE(3, 6, 4)
36 if (flg != PETSC_TRUE) {
38 "*** ERROR -my_file (MESH FILE NEEDED)");
54 std::cout <<
"<<<< NODESETs >>>>>" << std::endl;
57 std::cout << *it << std::endl;
58 CHKERR it->printBcData(std::cout);
59 std::vector<char> bc_data;
60 CHKERR it->getBcData(bc_data);
65 if (strcmp(&bc_data[0],
"Displacement") == 0) {
67 CHKERR it->getBcDataStructure(mydata);
74 else if (strcmp(&bc_data[0],
"Force") == 0) {
76 CHKERR it->getBcDataStructure(mydata);
83 else if (strcmp(&bc_data[0],
"Velocity") == 0) {
85 CHKERR it->getBcDataStructure(mydata);
92 else if (strcmp(&bc_data[0],
"Acceleration") == 0) {
94 CHKERR it->getBcDataStructure(mydata);
101 else if (strcmp(&bc_data[0],
"Temperature") == 0) {
103 CHKERR it->getBcDataStructure(mydata);
110 SETERRQ(PETSC_COMM_SELF, 1,
"Error: Unrecognizable BC type");
113 std::cout <<
"<<<< SIDESETs >>>>>" << std::endl;
116 std::cout << *it << std::endl;
117 CHKERR it->printBcData(std::cout);
118 std::vector<char> bc_data;
119 CHKERR it->getBcData(bc_data);
124 if (strcmp(&bc_data[0],
"Pressure") == 0) {
126 CHKERR it->getBcDataStructure(mydata);
133 else if (strcmp(&bc_data[0],
"HeatFlux") == 0) {
135 CHKERR it->getBcDataStructure(mydata);
142 else if (strcmp(&bc_data[0],
"cfd_bc") == 0) {
144 CHKERR it->getBcDataStructure(mydata);
147 if (mydata.
data.type == 6) {
150 std::cout << std::endl <<
"Interface" << std::endl;
151 myfile << std::endl <<
"Interface" << std::endl;
156 else if (mydata.
data.type == 15) {
160 std::cout << std::endl <<
"Pressure Inlet" << std::endl;
161 myfile << std::endl <<
"Pressure Inlet" << std::endl;
166 else if (mydata.
data.type == 16) {
170 std::cout << std::endl <<
"Pressure Outlet" << std::endl;
171 myfile << std::endl <<
"Pressure Outlet" << std::endl;
178 SETERRQ(PETSC_COMM_SELF, 1,
"Error: Unrecognizable BC type");
184 std::cout <<
"<<<< BLOCKSETs >>>>>" << std::endl;
187 std::cout << std::endl << *it << std::endl;
190 CHKERR it->printName(std::cout);
191 CHKERR it->printName(myfile);
194 std::vector<double> attributes;
195 CHKERR it->getAttributes(attributes);
196 CHKERR it->printAttributes(std::cout);
197 CHKERR it->printAttributes(myfile);
214 std::cout << std::endl << *it << std::endl;
217 std::string name = it->getName();
220 if (name.compare(0, 20,
"MAT_ELASTIC_TRANSISO") == 0) {
222 CHKERR it->getAttributeDataStructure(mydata);
226 }
else if (name.compare(0, 11,
"MAT_ELASTIC") == 0) {
228 CHKERR it->getAttributeDataStructure(mydata);
232 }
else if (name.compare(0, 10,
"MAT_INTERF") == 0) {
234 CHKERR it->getAttributeDataStructure(mydata);
239 SETERRQ(PETSC_COMM_SELF, 1,
"Error: Unrecognizable Material type");
#define CATCH_ERRORS
Catch errors.
#define CHKERR
Inline error check.
#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.
implementation of Data Operators for Forces and Sources
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition of the acceleration bc data structure.
Definition of the cfd_bc data structure.
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.
Definition of the displacement bc data structure.
Definition of the force bc data structure.
Definition of the heat flux bc data structure.
Transverse Isotropic material data structure.
Elastic material data structure.
Linear interface data structure.
Interface for managing meshsets containing materials and boundary conditions.
Definition of the pressure bc data structure.
Definition of the temperature bc data structure.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition of the velocity bc data structure.