Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_MeshFactory.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef Panzer_STK_MeshFactory_hpp__
12#define Panzer_STK_MeshFactory_hpp__
13
14#include <Teuchos_RCP.hpp>
15#include <Teuchos_ParameterListAcceptorDefaultBase.hpp>
16
17#include <stk_util/parallel/Parallel.hpp>
18
20
21namespace panzer_stk {
22
23class STK_Interface;
24
28class STK_MeshFactory : public Teuchos::ParameterListAcceptorDefaultBase {
29public:
31
40 virtual Teuchos::RCP<STK_Interface> buildMesh(stk::ParallelMachine parallelMach) const = 0;
41
46 virtual Teuchos::RCP<STK_Interface> buildUncommitedMesh(stk::ParallelMachine parallelMach) const = 0;
47
50 virtual void completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const = 0;
51
55 static void parsePeriodicBCList(const Teuchos::RCP<Teuchos::ParameterList> & pl,
56 std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > & periodicBC,
57 bool & useBBoxSearch)
58 {
60 parser.setParameterList(pl);
61 periodicBC = parser.getMatchers();
62 useBBoxSearch = parser.useBoundingBoxSearch();
63 }
64
65 void enableRebalance(bool enable,const Teuchos::RCP<const Teuchos::ParameterList> & rebalanceList=Teuchos::null)
66 { enableRebalance_ = enable;
67 rebalanceList_ = rebalanceList; }
68
69 void rebalance(STK_Interface & mesh) const
70 {
71 if(rebalanceList_!=Teuchos::null) {
72 // loop over user specified partitioning lists
73 for(Teuchos::ParameterList::ConstIterator itr=rebalanceList_->begin();
74 itr!=rebalanceList_->end();++itr) {
75
76 const Teuchos::ParameterEntry & entry = rebalanceList_->entry(itr);
77 TEUCHOS_TEST_FOR_EXCEPTION(!entry.isList(),std::runtime_error,
78 "Rebalance list is incorrect:\n" << entry << "\nA Zoltan list formated with strings is expected.");
79
80 // partition according to the list
81 mesh.rebalance(Teuchos::getValue<Teuchos::ParameterList>(entry));
82
83 // rebuild mesh internals
85 }
86 }
87 else if(enableRebalance_) {
88 // do the default thing, once
89 Teuchos::ParameterList emptyList;
90 mesh.rebalance(emptyList);
91
92 // rebuild mesh internals
94 }
95 }
96
97 double getMeshCoord(const int nx, const double deltaX, const double x0) const {
98 double x = static_cast<double>(nx)*deltaX;
99 double modX = std::abs(x);
100 double modX0 = std::abs(x0);
101 double val = x+x0;
102 if ((x0*x < 0.0) && (std::abs(modX-modX0) < std::numeric_limits<double>::epsilon()*modX0)) val=0.0;
103 return (val);
104 }
105
106protected:
107 // vector of periodic boundary condition objects
108 std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > periodicBCVec_;
109 // flag indicating which periodic search algorithm to use (bounding box or direct search)
111
112 // for managing rebalance
114 Teuchos::RCP<const Teuchos::ParameterList> rebalanceList_;
115};
116
117}
118
119#endif
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getMatchers() const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl)
void rebalance(const Teuchos::ParameterList &params)
virtual Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const =0
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const =0
void rebalance(STK_Interface &mesh) const
static void parsePeriodicBCList(const Teuchos::RCP< Teuchos::ParameterList > &pl, std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &periodicBC, bool &useBBoxSearch)
double getMeshCoord(const int nx, const double deltaX, const double x0) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const =0
Teuchos::RCP< const Teuchos::ParameterList > rebalanceList_
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCVec_
void enableRebalance(bool enable, const Teuchos::RCP< const Teuchos::ParameterList > &rebalanceList=Teuchos::null)