35 mesh->initialize(parallelMach);
47 machRank_ = stk::parallel_machine_rank(parallelMach);
48 machSize_ = stk::parallel_machine_size(parallelMach);
83 static RCP<Teuchos::ParameterList> defaultParams;
86 if(defaultParams == Teuchos::null) {
87 defaultParams = rcp(
new Teuchos::ParameterList);
89 defaultParams->set<
int>(
"X Elements",2);
90 defaultParams->set<
int>(
"Y Elements",2);
104 typedef shards::Quadrilateral<4> QuadTopo;
105 const CellTopologyData * ctd = shards::getCellTopologyData<QuadTopo>();
106 const CellTopologyData * side_ctd = shards::CellTopology(ctd).getBaseCellTopologyData(1,0);
110 for(
int bx=0;bx<2;bx++) {
112 std::stringstream ebPostfix;
113 ebPostfix <<
"-" <<
"0_" << bx;
136 else TEUCHOS_ASSERT(
false);
143 int myXElems_start = (
machRank_==0 ? 0 : 2);
144 int myXElems_end = (
machRank_==0 ? 1 : 3);
145 int myYElems_start = 0;
146 int myYElems_end = 1;
154 double deltaX = (xf_-x0_)/
double(totalXElems);
155 double deltaY = (yf_-y0_)/
double(totalYElems);
157 std::vector<double> coord(2,0.0);
160 for(
int nx=myXElems_start;nx<myXElems_end+1;++nx) {
161 coord[0] = double(nx)*deltaX+x0_;
162 for(
int ny=myYElems_start;ny<myYElems_end+1;++ny) {
163 coord[1] = double(ny)*deltaY+y0_;
165 mesh.
addNode(ny*(totalXElems+1)+nx+1,coord);
169 std::stringstream blockName;
170 blockName <<
"eblock-" << xBlock <<
"_" << yBlock;
174 for(
int nx=myXElems_start;nx<myXElems_end;++nx) {
175 for(
int ny=myYElems_start;ny<myYElems_end;++ny) {
176 stk::mesh::EntityId gid = totalXElems*ny+nx+1;
177 std::vector<stk::mesh::EntityId> nodes(4);
178 nodes[0] = nx+1+ny*(totalXElems+1);
179 nodes[1] = nodes[0]+1;
180 nodes[2] = nodes[1]+(totalXElems+1);
181 nodes[3] = nodes[2]-1;
191 unsigned int minElements =
nXElems_/size;
192 unsigned int extra =
nXElems_ - minElements*size;
194 TEUCHOS_ASSERT(minElements>0);
200 nume = minElements+1;
201 start = rank*(minElements+1);
205 start = extra*(minElements+1)+(rank-extra)*minElements;
208 return std::make_pair(start+
nXElems_*xBlock,nume);
215 return std::make_pair(start,
nYElems_);
222 std::size_t totalXElems =
nXElems_*2;
223 std::size_t totalYElems =
nYElems_*1;
226 stk::mesh::Part * left = mesh.
getSideset(
"left");
227 stk::mesh::Part * right = mesh.
getSideset(
"right");
228 stk::mesh::Part * top = mesh.
getSideset(
"top");
229 stk::mesh::Part * bottom = mesh.
getSideset(
"bottom");
231 std::vector<stk::mesh::Entity> localElmts;
235 std::vector<stk::mesh::Entity>::const_iterator itr;
236 for(itr=localElmts.begin();itr!=localElmts.end();++itr) {
237 stk::mesh::Entity element = (*itr);
241 ny = (gid-1) / totalXElems;
242 nx = gid-ny*totalXElems-1;
247 if(nx+1==totalXElems) {
274 if(ny+1==totalYElems) {
~MultiBlockMeshFactory()
Destructor.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶mList)
From ParameterListAcceptor.
void initializeWithDefaults()
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
std::pair< int, int > determineXElemSizeAndStart(int xBlock, unsigned int size, unsigned int rank) const
void addSideSets(STK_Interface &mesh) const
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
From ParameterListAcceptor.
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void buildBlock(stk::ParallelMachine parallelMach, int xBlock, int yBlock, STK_Interface &mesh) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
MultiBlockMeshFactory()
Constructor.
std::pair< int, int > determineYElemSizeAndStart(int yBlock, unsigned int size, unsigned int rank) const
void buildLocalElementIDs()
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
bool isInitialized() const
Has initialize been called on this mesh object?
void buildSubcells()
force the mesh to build subcells: edges and faces
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addSideset(const std::string &name, const CellTopologyData *ctData)
void endModification(const bool find_and_set_shared_nodes_in_stk=true)
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
stk::mesh::Part * getSideset(const std::string &name) const
void rebalance(STK_Interface &mesh) const