v0.14.0
add_cubit_meshsets.cpp

Test and example setting cubit meshsets

/**
* @file add_cubit_meshsets.cpp
* @example add_cubit_meshsets.cpp
* @brief Test and example setting cubit meshsets
*
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::setLog("ATOM_TEST");
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
// Create MoFEM (Joseph) database
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
MeshsetsManager *meshsets_manager_ptr;
CHKERR m_field.getInterface(meshsets_manager_ptr);
MOFEM_LOG_CHANNEL("ATOM_TEST")
MOFEM_LOG_ATTRIBUTES("ATOM_TEST",
MOFEM_LOG_TAG("ATOM_TEST", "atom test");
MOFEM_LOG("ATOM_TEST", Sev::verbose) << "<<<< SIDESETs >>>>>";
bool add_block_is_there = false;
CHKERR meshsets_manager_ptr->addMeshset(SIDESET, 1002);
{
strncpy(mybc.data.name, "Pressure", 8);
mybc.data.flag1 = 0;
mybc.data.flag2 = 0;
mybc.data.value1 = 1;
CHKERR meshsets_manager_ptr->setBcData(SIDESET, 1002, mybc);
}
if (it->getMeshsetId() != 1002)
continue;
add_block_is_there = true;
CHKERR it->getBcDataStructure(mydata);
// Print data
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if (!add_block_is_there)
SETERRQ(PETSC_COMM_WORLD, MOFEM_OPERATION_UNSUCCESSFUL,
"no added block set");
MOFEM_LOG("ATOM_TEST", Sev::inform) << "<<<< BLOCKSETs >>>>>";
add_block_is_there = false;
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1000, "ADD_BLOCK_SET");
std::vector<double> attr(3);
attr[0] = 0;
attr[1] = 1;
attr[2] = 2;
CHKERR meshsets_manager_ptr->setAtributes(BLOCKSET, 1000, attr);
// Get block name
std::string name = it->getName();
if (name.compare(0, 13, "ADD_BLOCK_SET") == 0) {
add_block_is_there = true;
std::vector<double> attributes;
it->getAttributes(attributes);
if (attributes.size() != 3) {
SETERRQ1(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"should be 3 attributes but is %d", attributes.size());
}
if (attributes[0] != 0 || attributes[1] != 1 || attributes[2] != 2) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"wrong values of attributes");
}
}
}
if (!add_block_is_there) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_OPERATION_UNSUCCESSFUL,
"no added block set");
}
add_block_is_there = false;
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1001, "MAT_ELASTIC");
{
Mat_Elastic mydata;
mydata.data.Young = 1;
mydata.data.Poisson = 0.25;
CHKERR meshsets_manager_ptr->setAtributesByDataStructure(BLOCKSET, 1001,
mydata);
}
if (it->getMeshsetId() != 1001)
continue;
// Get block name
std::string name = it->getName();
if (name.compare(0, 13, "MAT_ELASTIC") == 0 &&
(it->getBcType() & CubitBCType(MAT_ELASTICSET)).any()) {
add_block_is_there = true;
Mat_Elastic mydata;
CHKERR it->getAttributeDataStructure(mydata);
// Print data
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
if (mydata.data.Young != 1 || mydata.data.Poisson != 0.25) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"wrong values of attributes");
}
}
}
if (!add_block_is_there) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_OPERATION_UNSUCCESSFUL,
"no added block set");
}
MOFEM_LOG("ATOM_TEST", Sev::inform) << "<<<< NODESET >>>>>";
CHKERR meshsets_manager_ptr->addMeshset(NODESET, 1010);
std::memcpy(disp_bc.data.name, "Displacement", 12);
disp_bc.data.flag1 = 1;
disp_bc.data.flag2 = 1;
disp_bc.data.flag3 = 1;
disp_bc.data.flag4 = 0;
disp_bc.data.flag5 = 0;
disp_bc.data.flag6 = 0;
disp_bc.data.value1 = 0;
disp_bc.data.value2 = 0;
disp_bc.data.value3 = 0;
disp_bc.data.value4 = 0;
disp_bc.data.value5 = 0;
disp_bc.data.value6 = 0;
CHKERR meshsets_manager_ptr->setBcData(NODESET, 1010, disp_bc);
m_field, NODESET | DISPLACEMENTSET, it)) {
CHKERR it->getBcDataStructure(disp_data);
MOFEM_LOG("ATOM_TEST", Sev::inform) << disp_data;
}
MOFEM_LOG("ATOM_TEST", Sev::inform)
<< "<<<< ADD BLOCKSETs FROM CONFIG FILE >>>>>";
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1002, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1003, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1004, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1005, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1006, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1007, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1008, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1009, "ADD_BLOCK_SET");
CHKERR meshsets_manager_ptr->addMeshset(BLOCKSET, 1010, "foo");
CHKERR meshsets_manager_ptr->setMeshsetFromFile(
/*"add_cubit_meshsets.in"*/);
// List all meshsets
MOFEM_LOG("ATOM_TEST", Sev::inform) << "Iterate blocksets";
bool mat_elastic_trans_is_found = true;
for (_IT_CUBITMESHSETS_FOR_LOOP_(m_field, it)) {
MOFEM_LOG("ATOM_TEST", Sev::inform) << *it;
if ((it->getBcType() & CubitBCType(BLOCKSET)).any()) {
std::vector<double> attributes;
it->getAttributes(attributes);
std::ostringstream ss;
ss << "Attr: ";
for (unsigned int ii = 0; ii != attributes.size(); ii++) {
ss << attributes[ii] << " ";
}
MOFEM_LOG("ATOM_TEST", Sev::inform) << ss.str();
std::string block_name = it->getName();
if (block_name.compare(0, block_name.size(), "MAT_ELASTIC_TRANS_ISO") ==
0) {
MOFEM_LOG("ATOM_TEST", Sev::inform) << "Mat Trans Iso block ";
mat_elastic_trans_is_found = true;
CHKERR it->getAttributeDataStructure(mydata);
// Print data
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
if (mydata.data.Youngp != 1)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value");
if (mydata.data.Youngz != 2)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value");
if (mydata.data.Poissonp != 3)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value");
if (mydata.data.Poissonpz != 4)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value");
if (mydata.data.Shearzp != 5)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Wrong value");
}
}
if(!mat_elastic_trans_is_found)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "Block not found");
if ((it->getBcType() & CubitBCType(MAT_ELASTICSET)).any()) {
Mat_Elastic mydata;
CHKERR it->getAttributeDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << "Mat elastic found ";
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(MAT_THERMALSET)).any()) {
Mat_Thermal mydata;
CHKERR it->getAttributeDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << "Mat thermal found ";
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(DISPLACEMENTSET)).any()) {
CHKERR it->getBcDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(FORCESET)).any()) {
CHKERR it->getBcDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(PRESSURESET)).any()) {
CHKERR it->getBcDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(TEMPERATURESET)).any()) {
CHKERR it->getBcDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
if ((it->getBcType() & CubitBCType(HEATFLUXSET)).any()) {
CHKERR it->getBcDataStructure(mydata);
MOFEM_LOG("ATOM_TEST", Sev::inform) << mydata;
}
}
}
}
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
SIDESET
@ SIDESET
Definition: definitions.h:160
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
MoFEM::TemperatureCubitBcData
Definition of the temperature bc data structure.
Definition: BCData.hpp:306
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
MoFEM::Mat_Elastic
Elastic material data structure.
Definition: MaterialBlocks.hpp:139
_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Definition: MeshsetsManager.hpp:49
MOFEM_LOG_ATTRIBUTES
#define MOFEM_LOG_ATTRIBUTES(channel, bit)
Add attributes to channel.
Definition: LogManager.hpp:296
MoFEM.hpp
MoFEM::ForceCubitBcData
Definition of the force bc data structure.
Definition: BCData.hpp:139
MoFEM::DisplacementCubitBcData
Definition of the displacement bc data structure.
Definition: BCData.hpp:76
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
MoFEM::PressureCubitBcData
Definition of the pressure bc data structure.
Definition: BCData.hpp:379
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
NODESET
@ NODESET
Definition: definitions.h:159
MoFEM::LogManager::BitLineID
@ BitLineID
Definition: LogManager.hpp:48
MoFEM::LogManager::BitScope
@ BitScope
Definition: LogManager.hpp:49
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
main
int main(int argc, char *argv[])
Definition: add_cubit_meshsets.cpp:14
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
FORCESET
@ FORCESET
Definition: definitions.h:164
MAT_THERMALSET
@ MAT_THERMALSET
block name is "MAT_THERMAL"
Definition: definitions.h:174
MoFEM::MeshsetsManager::setMeshsetFromFile
MoFEMErrorCode setMeshsetFromFile(const string file_name, const bool clean_file_options=true)
add blocksets reading config file
Definition: MeshsetsManager.cpp:791
MoFEM::Mat_Elastic_TransIso::data
_data_ data
Definition: MaterialBlocks.hpp:382
MoFEM::Mat_Elastic_TransIso
Transverse Isotropic material data structure.
Definition: MaterialBlocks.hpp:371
MoFEM::MeshsetsManager::addMeshset
MoFEMErrorCode addMeshset(const CubitBCType cubit_bc_type, const int ms_id, const std::string name="")
add cubit meshset
Definition: MeshsetsManager.cpp:385
MAT_ELASTICSET
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:172
DISPLACEMENTSET
@ DISPLACEMENTSET
Definition: definitions.h:163
MoFEM::MeshsetsManager::setAtributesByDataStructure
MoFEMErrorCode setAtributesByDataStructure(const CubitBCType cubit_bc_type, const int ms_id, const GenericAttributeData &data, const std::string name="")
set (material) data structure to cubit meshset
Definition: MeshsetsManager.cpp:499
MoFEM::LogManager::getStrmSelf
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
Definition: LogManager.cpp:340
MoFEM::Mat_Elastic::data
_data_ data
Definition: MaterialBlocks.hpp:155
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::DisplacementCubitBcData::data
_data_ data
Definition: BCData.hpp:99
MoFEM::MeshsetsManager::setAtributes
MoFEMErrorCode setAtributes(const CubitBCType cubit_bc_type, const int ms_id, const std::vector< double > &attributes, const std::string name="")
set attributes to cubit meshset
Definition: MeshsetsManager.cpp:459
MoFEM::PressureCubitBcData::data
_data_ data
Definition: BCData.hpp:387
MoFEM::CoreTmp< 0 >::Initialize
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
MoFEM::Types::CubitBCType
std::bitset< 32 > CubitBCType
Definition: Types.hpp:52
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#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.
Definition: MeshsetsManager.hpp:71
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
MoFEM::Mat_Thermal
Thermal material data structure.
Definition: MaterialBlocks.hpp:201
help
static char help[]
Definition: add_cubit_meshsets.cpp:12
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
TEMPERATURESET
@ TEMPERATURESET
Definition: definitions.h:168
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::MeshsetsManager::setBcData
MoFEMErrorCode setBcData(const CubitBCType cubit_bc_type, const int ms_id, const GenericCubitBcData &data)
set boundary data structure to meshset
Definition: MeshsetsManager.cpp:530
MoFEM::HeatFluxCubitBcData
Definition of the heat flux bc data structure.
Definition: BCData.hpp:427
_IT_CUBITMESHSETS_FOR_LOOP_
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
Definition: MeshsetsManager.hpp:34
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
HEATFLUXSET
@ HEATFLUXSET
Definition: definitions.h:169
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389