13 #include <PanzerAdaptersSTK_config.hpp>
18 namespace panzer_stk {
33 PANZER_FUNC_TIME_MONITOR(
"panzer::CubeTetMeshFactory::buildMesh()");
39 mesh->initialize(parallelMach);
49 PANZER_FUNC_TIME_MONITOR(
"panzer::CubeTetMeshFactory::buildUncomittedMesh()");
53 machRank_ = stk::parallel_machine_rank(parallelMach);
54 machSize_ = stk::parallel_machine_size(parallelMach);
66 const int maxFactor = 50;
69 int factors[maxFactor];
70 for (
int jj = 0; jj < maxFactor; jj++) factors[jj] = 0;
71 for (
int jj = 2; jj < maxFactor; jj++) {
74 int temp = ProcTemp/jj;
75 if (temp*jj == ProcTemp) {
85 for (
int jj = maxFactor-1; jj > 0; jj--) {
86 while (factors[jj] != 0) {
102 "Cannot build CubeTetMeshFactory, the product of \"X Procs\", \"Y Procs\", and \"Z Procs\""
103 " must equal the number of processors.");
117 PANZER_FUNC_TIME_MONITOR(
"panzer::CubeTetMeshFactory::completeMeshConstruction()");
159 x0_ = paramList->
get<
double>(
"X0");
160 y0_ = paramList->
get<
double>(
"Y0");
161 z0_ = paramList->
get<
double>(
"Z0");
163 xf_ = paramList->
get<
double>(
"Xf");
164 yf_ = paramList->
get<
double>(
"Yf");
165 zf_ = paramList->
get<
double>(
"Zf");
192 if(defaultParams == Teuchos::null) {
195 defaultParams->
set<
double>(
"X0",0.0);
196 defaultParams->
set<
double>(
"Y0",0.0);
197 defaultParams->
set<
double>(
"Z0",0.0);
199 defaultParams->
set<
double>(
"Xf",1.0);
200 defaultParams->
set<
double>(
"Yf",1.0);
201 defaultParams->
set<
double>(
"Zf",1.0);
203 defaultParams->
set<
int>(
"X Blocks",1);
204 defaultParams->
set<
int>(
"Y Blocks",1);
205 defaultParams->
set<
int>(
"Z Blocks",1);
207 defaultParams->
set<
int>(
"X Procs",-1);
208 defaultParams->
set<
int>(
"Y Procs",1);
209 defaultParams->
set<
int>(
"Z Procs",1);
211 defaultParams->
set<
int>(
"X Elements",5);
212 defaultParams->
set<
int>(
"Y Elements",5);
213 defaultParams->
set<
int>(
"Z Elements",5);
216 defaultParams->
set<
bool>(
"Create Edge Blocks",
false,
"Create edge blocks in the mesh");
217 defaultParams->
set<
bool>(
"Create Face Blocks",
false,
"Create face blocks in the mesh");
220 bcs.
set<
int>(
"Count",0);
223 return defaultParams;
245 typedef shards::Tetrahedron<4> TetTopo;
246 const CellTopologyData * ctd = shards::getCellTopologyData<TetTopo>();
247 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(2,0);
248 const CellTopologyData * edge_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
249 const CellTopologyData * face_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(2,0);
257 std::stringstream ebPostfix;
258 ebPostfix <<
"-" << bx <<
"_" << by <<
"_" << bz;
291 for(
int xBlock=0;xBlock<
xBlocks_;xBlock++) {
292 for(
int yBlock=0;yBlock<
yBlocks_;yBlock++) {
293 for(
int zBlock=0;zBlock<
zBlocks_;zBlock++) {
294 buildBlock(parallelMach,xBlock,yBlock,zBlock,mesh);
308 int myXElems_start = sizeAndStartX.first;
309 int myXElems_end = myXElems_start+sizeAndStartX.second;
310 int myYElems_start = sizeAndStartY.first;
311 int myYElems_end = myYElems_start+sizeAndStartY.second;
312 int myZElems_start = sizeAndStartZ.first;
313 int myZElems_end = myZElems_start+sizeAndStartZ.second;
319 double deltaX = (
xf_-
x0_)/
double(totalXElems);
320 double deltaY = (
yf_-
y0_)/
double(totalYElems);
321 double deltaZ = (
zf_-
z0_)/
double(totalZElems);
323 std::vector<double> coord(3,0.0);
326 for(
int nx=myXElems_start;nx<myXElems_end+1;++nx) {
328 for(
int ny=myYElems_start;ny<myYElems_end+1;++ny) {
330 for(
int nz=myZElems_start;nz<myZElems_end+1;++nz) {
333 mesh.
addNode(nz*(totalYElems+1)*(totalXElems+1)+ny*(totalXElems+1)+nx+1,coord);
338 std::stringstream blockName;
339 blockName <<
"eblock-" << xBlock <<
"_" << yBlock <<
"_" << zBlock;
343 for(
int nx=myXElems_start;nx<myXElems_end;++nx) {
344 for(
int ny=myYElems_start;ny<myYElems_end;++ny) {
345 for(
int nz=myZElems_start;nz<myZElems_end;++nz) {
347 std::vector<stk::mesh::EntityId> nodes(8);
348 nodes[0] = nx+1+ny*(totalXElems+1) +nz*(totalYElems+1)*(totalXElems+1);
349 nodes[1] = nodes[0]+1;
350 nodes[2] = nodes[1]+(totalXElems+1);
351 nodes[3] = nodes[2]-1;
352 nodes[4] = nodes[0]+(totalYElems+1)*(totalXElems+1);
353 nodes[5] = nodes[1]+(totalYElems+1)*(totalXElems+1);
354 nodes[6] = nodes[2]+(totalYElems+1)*(totalXElems+1);
355 nodes[7] = nodes[3]+(totalYElems+1)*(totalXElems+1);
357 buildTetsOnHex(Teuchos::tuple(totalXElems,totalYElems,totalZElems),
358 Teuchos::tuple(nx,ny,nz),
367 stk::mesh::Part * block,
368 const std::vector<stk::mesh::EntityId> & h_nodes,
375 int totalXElems = meshDesc[0];
int totalYElems = meshDesc[1];
int totalZElems = meshDesc[2];
376 int nx = element[0];
int ny = element[1];
int nz = element[2];
378 stk::mesh::EntityId hex_id = totalXElems*totalYElems*nz+totalXElems*ny+nx+1;
379 stk::mesh::EntityId gid_0 = 12*(hex_id-1)+1;
380 std::vector<stk::mesh::EntityId> nodes(4);
383 stk::mesh::EntityId centroid = 0;
385 stk::mesh::EntityId largestNode = (totalXElems+1)*(totalYElems+1)*(totalZElems+1);
386 centroid = hex_id+largestNode;
389 std::vector<double> coord(3,0.0);
390 for(std::size_t i=0;i<h_nodes.size();i++) {
392 coord[0] += node_coord[0];
393 coord[1] += node_coord[1];
394 coord[2] += node_coord[2];
404 int idSet[][3] = { { 0, 1, 2},
417 for(
int i=0;i<12;i++) {
418 nodes[0] = h_nodes[idSet[i][0]];
419 nodes[1] = h_nodes[idSet[i][1]];
420 nodes[2] = h_nodes[idSet[i][2]];
431 unsigned int minElements =
nXElems_/size;
432 unsigned int extra =
nXElems_ - minElements*size;
440 nume = minElements+1;
441 start = xProcLoc*(minElements+1);
445 start = extra*(minElements+1)+(xProcLoc-extra)*minElements;
448 return std::make_pair(start+
nXElems_*xBlock,nume);
457 unsigned int minElements =
nYElems_/size;
458 unsigned int extra =
nYElems_ - minElements*size;
466 nume = minElements+1;
467 start = yProcLoc*(minElements+1);
471 start = extra*(minElements+1)+(yProcLoc-extra)*minElements;
474 return std::make_pair(start+
nYElems_*yBlock,nume);
482 unsigned int minElements =
nZElems_/size;
483 unsigned int extra =
nZElems_ - minElements*size;
491 nume = minElements+1;
492 start = zProcLoc*(minElements+1);
496 start = extra*(minElements+1)+(zProcLoc-extra)*minElements;
499 return std::make_pair(start+
nZElems_*zBlock,nume);
505 const stk::mesh::EntityRank side_rank = mesh.
getSideRank();
512 stk::mesh::Part * left = mesh.
getSideset(
"left");
513 stk::mesh::Part * right = mesh.
getSideset(
"right");
514 stk::mesh::Part * top = mesh.
getSideset(
"top");
515 stk::mesh::Part * bottom = mesh.
getSideset(
"bottom");
516 stk::mesh::Part * front = mesh.
getSideset(
"front");
517 stk::mesh::Part * back = mesh.
getSideset(
"back");
519 std::vector<stk::mesh::Entity> localElmts;
525 std::vector<stk::mesh::Entity>::const_iterator itr;
526 for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
527 stk::mesh::Entity element = (*itr);
531 stk::mesh::EntityId h_gid = (gid-1)/12+1;
532 stk::mesh::EntityId t_offset = gid - (12*(h_gid-1)+1);
534 std::size_t nx,ny,nz;
535 nz = (h_gid-1) / (totalXElems*totalYElems);
536 h_gid = (h_gid-1)-nz*(totalXElems*totalYElems);
537 ny = h_gid / totalXElems;
538 nx = h_gid-ny*totalXElems;
540 if(nz==0 && (t_offset==0 || t_offset==1)) {
547 if(nz+1==totalZElems && (t_offset==10 || t_offset==11)) {
555 if(ny==0 && (t_offset==2 || t_offset==3)) {
562 if(ny+1==totalYElems && (t_offset==8 || t_offset==9)) {
570 if(nx==0 && (t_offset==4 || t_offset==5)) {
577 if(nx+1==totalXElems && (t_offset==6 || t_offset==7)) {
594 stk::mesh::Part * origin = mesh.
getNodeset(
"origin");
600 stk::mesh::Entity node = bulkData->get_entity(mesh.
getNodeRank(),1);
616 stk::mesh::Selector owned_block = metaData->locally_owned_part();
618 std::vector<stk::mesh::Entity> edges;
619 bulkData->get_entities(mesh.
getEdgeRank(), owned_block, edges);
632 stk::mesh::Selector owned_block = metaData->locally_owned_part();
634 std::vector<stk::mesh::Entity> faces;
635 bulkData->get_entities(mesh.
getFaceRank(), owned_block, faces);
642 std::size_t i=0,j=0,k=0;
648 return Teuchos::tuple(i,j,k);
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
std::pair< int, int > determineYElemSizeAndStart(int yBlock, unsigned int size, unsigned int rank) const
void addNodeset(const std::string &name)
std::string edgeBlockName_
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
basic_FancyOStream & setShowProcRank(const bool showProcRank)
T & get(const std::string &name, T def_value)
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
std::string faceBlockName_
void addFaceBlocks(STK_Interface &mesh) const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
std::pair< int, int > determineZElemSizeAndStart(int zBlock, unsigned int size, unsigned int rank) const
stk::mesh::Part * getNodeset(const std::string &name) const
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
stk::mesh::Part * getSideset(const std::string &name) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void addSideSets(STK_Interface &mesh) const
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
stk::mesh::EntityRank getSideRank() const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
void rebalance(STK_Interface &mesh) const
void buildBlock(stk::ParallelMachine machRank, int xBlock, int yBlock, int zBlock, STK_Interface &mesh) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
void setMyParamList(const RCP< ParameterList > ¶mList)
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
static const std::string edgeBlockString
virtual ~CubeTetMeshFactory()
Destructor.
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶mList)
From ParameterListAcceptor.
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
stk::mesh::EntityRank getFaceRank() const
void addEdgeBlocks(STK_Interface &mesh) const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
bool isInitialized() const
Has initialize been called on this mesh object?
Teuchos::Tuple< std::size_t, 3 > procRankToProcTuple(std::size_t procRank) const
what is the 3D tuple describe this processor distribution
void initializeWithDefaults()
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
static const std::string faceBlockString
double getMeshCoord(const int nx, const double deltaX, const double x0) const
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
void buildTetsOnHex(const Teuchos::Tuple< int, 3 > &meshDesc, const Teuchos::Tuple< int, 3 > &element, stk::mesh::Part *block, const std::vector< stk::mesh::EntityId > &h_nodes, STK_Interface &mesh) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
stk::mesh::EntityRank getEdgeRank() const
void addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
void buildLocalElementIDs()
void addFaceBlock(const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
#define TEUCHOS_ASSERT(assertion_test)
stk::mesh::EntityRank getNodeRank() const
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
CubeTetMeshFactory()
Constructor.
Teuchos::Tuple< std::size_t, 3 > procTuple_
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
void addNodeSets(STK_Interface &mesh) const