11 #ifndef __Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__
12 #define __Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__
15 #include "Panzer_STK_ScatterFields.hpp"
16 #include "Panzer_STK_ScatterVectorFields.hpp"
17 #include "Panzer_PointValues_Evaluator.hpp"
18 #include "Panzer_BasisValues_Evaluator.hpp"
22 #include <unordered_set>
24 namespace panzer_stk {
30 Response_STKDummy(
const std::string & rn)
32 virtual void scatterResponse() {}
33 virtual void initializeResponse() {}
36 Response_STKDummy(
const Response_STKDummy &);
40 template <
typename EvalT>
44 return Teuchos::rcp(
new Response_STKDummy(responseName));
47 template <
typename EvalT>
57 typedef std::pair<std::string,RCP<const panzer::PureBasis> > StrConstPureBasisPair;
60 std::unordered_set<std::string> scaledFieldsHash = scaledFieldsHash_;
62 std::map<std::string,RCP<const panzer::PureBasis> > bases;
63 std::map<std::string,std::vector<std::string> > basisBucket;
65 const std::map<std::string,RCP<panzer::PureBasis> > & nc_bases = physicsBlock.
getBases();
66 bases.insert(nc_bases.begin(),nc_bases.end());
69 std::vector<StrConstPureBasisPair> allFields;
73 if(!addCoordinateFields_ && addSolutionFields_) {
79 std::vector<std::string> removedFields;
80 const std::vector<std::vector<std::string> > & coord_fields = physicsBlock.
getCoordinateDOFs();
81 for(std::size_t c=0;c<coord_fields.size();c++)
82 for(std::size_t d=0;d<coord_fields[c].size();d++)
83 removedFields.push_back(coord_fields[c][d]);
86 deleteRemovedFields(removedFields,allFields);
88 else if(addCoordinateFields_ && !addSolutionFields_) {
90 const std::vector<std::vector<std::string> > & coord_fields = physicsBlock.
getCoordinateDOFs();
93 for(std::size_t c=0;c<coord_fields.size();c++) {
94 for(std::size_t d=0;d<coord_fields[c].size();d++) {
99 allFields.push_back(std::make_pair(coord_fields[c][d],basis));
103 else if(addSolutionFields_)
107 if(addSolutionFields_)
111 for(std::size_t i=0;i<additionalFields_.size();i++)
112 bases[additionalFields_[i].second->name()] = additionalFields_[i].second;
114 allFields.insert(allFields.end(),additionalFields_.begin(),additionalFields_.end());
116 deleteRemovedFields(removedFields_,allFields);
118 bucketByBasisType(allFields,basisBucket);
124 itr!=bases.end();++itr) {
126 if(itr->second->isVectorBasis()) {
130 Kokkos::DynRankView<double,PHX::Device> centroid;
136 this->
template registerEvaluator<EvalT>(fm, evaluator);
145 for(std::map<std::string,std::vector<std::string> >::const_iterator itr=basisBucket.begin();
146 itr!=basisBucket.end();++itr) {
148 std::string basisName = itr->first;
149 const std::vector<std::string> & fields = itr->second;
151 std::map<std::string,Teuchos::RCP<const panzer::PureBasis> >::const_iterator found = bases.find(basisName);
153 "Could not find basis \""+basisName+
"\"!");
157 std::vector<double> scalars(fields.size(),1.0);
158 for(std::size_t f=0;f<fields.size();f++) {
159 std::unordered_map<std::string,double>::const_iterator f2s_itr = fieldToScalar_.find(fields[f]);
163 if(f2s_itr!=fieldToScalar_.end()) {
164 scalars[f] = f2s_itr->second;
165 scaledFieldsHash.erase(fields[f]);
173 std::string fields_concat =
"";
174 for(std::size_t f=0;f<fields.size();f++) {
175 fields_concat += fields[f];
180 mesh_, basis, fields,scalars));
183 this->
template registerEvaluator<EvalT>(fm, eval);
184 fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
193 this->
template registerEvaluator<EvalT>(fm, evaluator);
197 std::string fields_concat =
"";
198 for(std::size_t f=0;f<fields.size();f++) {
200 p.
set(
"Name",fields[f]);
201 p.
set(
"Basis",basis);
206 this->
template registerEvaluator<EvalT>(fm, evaluator);
208 fields_concat += fields[f];
215 mesh_,centroidRule,fields,scalars));
217 this->
template registerEvaluator<EvalT>(fm, evaluator);
218 fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]);
228 this->
template registerEvaluator<EvalT>(fm, evaluator);
232 std::string fields_concat =
"";
233 for(std::size_t f=0;f<fields.size();f++) {
235 p.
set(
"Name",fields[f]);
236 p.
set(
"Basis",basis);
241 this->
template registerEvaluator<EvalT>(fm, evaluator);
243 fields_concat += fields[f];
250 mesh_,centroidRule,fields,scalars));
252 this->
template registerEvaluator<EvalT>(fm, evaluator);
253 fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]);
262 for(std::unordered_set<std::string>::const_iterator itr=scaledFieldsHash.begin();
263 itr!=scaledFieldsHash.end();itr++) {
264 out <<
"WARNING: STK Solution Writer did not scale the field \"" << *itr <<
"\" "
265 <<
"because it was not written." << std::endl;
269 template <
typename EvalT>
272 std::map<std::string,std::vector<std::string> > & basisBucket)
275 for(std::size_t i=0;i<providedDofs.size();i++) {
276 std::string fieldName = providedDofs[i].first;
279 basisBucket[basis->
name()].push_back(fieldName);
283 template <
typename EvalT>
287 Kokkos::DynRankView<double,PHX::Device> & centroid)
const
290 using Teuchos::rcp_dynamic_cast;
292 centroid = Kokkos::DynRankView<double,PHX::Device>(
"centroid",1,baseDimension);
293 auto l_centroid = centroid;
297 itr!=bases.end();++itr) {
302 Kokkos::DynRankView<double,PHX::Device> coords(
"coords",intrepidBasis->getCardinality(),
303 intrepidBasis->getBaseCellTopology().getDimension());
304 intrepidBasis->getDofCoords(coords);
308 Kokkos::parallel_for(coords.extent_int(0), KOKKOS_LAMBDA (
int i) {
309 for(
int d=0;d<coords.extent_int(1);d++)
310 l_centroid(0,d) += coords(i,d)/coords.extent(0);
321 template <
typename EvalT>
325 fieldToScalar_[fieldName] = fieldScalar;
326 scaledFieldsHash_.insert(fieldName);
329 template <
typename EvalT>
333 if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>())
339 template <
typename EvalT>
343 additionalFields_.push_back(std::make_pair(fieldName,basis));
346 template <
typename EvalT>
355 fields.erase(std::remove_if(fields.begin(),fields.end(),functor),fields.end());
std::string name() const
A unique key that is the combination of the basis type and basis order.
void deleteRemovedFields(const std::vector< std::string > &removedFields, std::vector< std::pair< std::string, Teuchos::RCP< const panzer::PureBasis > > > &fields) const
Delete from the argument all the fields that are in the removedFields array.
RCP< const T > getConst() const
int baseCellDimension() const
Dimension of the base cell. NOT the dimension of the local side, even if the side() method returns tr...
Object that contains information on the physics and discretization of a block of elements with the SA...
const std::map< std::string, Teuchos::RCP< panzer::PureBasis > > & getBases() const
Returns the unique set of bases, key is the unique panzer::PureBasis::name() of the basis...
const panzer::CellData & cellData() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual Teuchos::RCP< panzer::ResponseBase > buildResponseObject(const std::string &responseName) const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Interpolates basis DOF values to IP DOF values.
virtual bool typeSupported() const
void addAdditionalField(const std::string &fieldName, const Teuchos::RCP< const panzer::PureBasis > &basis)
Interpolates basis DOF values to IP DOF Curl values.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual Teuchos::RCP< const panzer::PureBasis > lookupBasis(const std::string &fieldName) const =0
Get the basis associated with a particular field.
const std::vector< std::vector< std::string > > & getCoordinateDOFs() const
void scaleField(const std::string &fieldName, double fieldScalar)
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
static void bucketByBasisType(const std::vector< std::pair< std::string, Teuchos::RCP< const panzer::PureBasis > > > &providedDofs, std::map< std::string, std::vector< std::string > > &basisBucket)
const std::vector< StrPureBasisPair > & getTangentFields() const
Returns list of tangent fields from DOFs and tangent param names.
Teuchos::RCP< const FieldLibraryBase > getFieldLibraryBase() const
virtual void buildAndRegisterEvaluators(const std::string &responseName, PHX::FieldManager< panzer::Traits > &fm, const panzer::PhysicsBlock &physicsBlock, const Teuchos::ParameterList &user_data) const
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
Interpolates basis DOF values to IP DOF values.
std::vector< std::string > removedFields_
const std::vector< StrPureBasisPair > & getProvidedDOFs() const
void computeReferenceCentroid(const std::map< std::string, Teuchos::RCP< const panzer::PureBasis > > &bases, int baseDimension, Kokkos::DynRankView< double, PHX::Device > ¢roid) const