1 #ifndef __Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__ 
    2 #define __Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__ 
    5 #include "Panzer_STK_ScatterFields.hpp" 
    6 #include "Panzer_STK_ScatterVectorFields.hpp" 
    7 #include "Panzer_PointValues_Evaluator.hpp" 
    8 #include "Panzer_BasisValues_Evaluator.hpp" 
   12 #include <unordered_set> 
   14 namespace panzer_stk {
 
   20      Response_STKDummy(
const std::string & rn)
 
   22      virtual void scatterResponse() {}
 
   23      virtual void initializeResponse() {}
 
   26      Response_STKDummy(
const Response_STKDummy &);
 
   30 template <
typename EvalT> 
 
   34    return Teuchos::rcp(
new Response_STKDummy(responseName));
 
   37 template <
typename EvalT> 
 
   47   typedef std::pair<std::string,RCP<const panzer::PureBasis> > StrConstPureBasisPair;
 
   50   std::unordered_set<std::string> scaledFieldsHash = scaledFieldsHash_;
 
   52   std::map<std::string,RCP<const panzer::PureBasis> > bases;
 
   53   std::map<std::string,std::vector<std::string> > basisBucket;
 
   55     const std::map<std::string,RCP<panzer::PureBasis> > & nc_bases = physicsBlock.
getBases();
 
   56     bases.insert(nc_bases.begin(),nc_bases.end());
 
   59   std::vector<StrConstPureBasisPair> allFields;
 
   63   if(!addCoordinateFields_ && addSolutionFields_) {
 
   69     std::vector<std::string> removedFields;
 
   70     const std::vector<std::vector<std::string> > & coord_fields = physicsBlock.
getCoordinateDOFs();
 
   71     for(std::size_t c=0;c<coord_fields.size();c++)
 
   72       for(std::size_t d=0;d<coord_fields[c].size();d++)
 
   73         removedFields.push_back(coord_fields[c][d]);
 
   76     deleteRemovedFields(removedFields,allFields); 
 
   78   else if(addCoordinateFields_ && !addSolutionFields_) {
 
   80     const std::vector<std::vector<std::string> > & coord_fields = physicsBlock.
getCoordinateDOFs();
 
   83     for(std::size_t c=0;c<coord_fields.size();c++) {
 
   84       for(std::size_t d=0;d<coord_fields[c].size();d++) {
 
   89         allFields.push_back(std::make_pair(coord_fields[c][d],basis));
 
   93   else if(addSolutionFields_)
 
   97   if(addSolutionFields_)
 
  101   for(std::size_t i=0;i<additionalFields_.size();i++)
 
  102     bases[additionalFields_[i].second->name()] = additionalFields_[i].second;
 
  104   allFields.insert(allFields.end(),additionalFields_.begin(),additionalFields_.end());
 
  106   deleteRemovedFields(removedFields_,allFields);
 
  108   bucketByBasisType(allFields,basisBucket);
 
  114       itr!=bases.end();++itr) {
 
  116     if(itr->second->isVectorBasis()) {
 
  120       Kokkos::DynRankView<double,PHX::Device> centroid;
 
  126       this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  135   for(std::map<std::string,std::vector<std::string> >::const_iterator itr=basisBucket.begin();
 
  136       itr!=basisBucket.end();++itr) {
 
  138     std::string basisName = itr->first;
 
  139     const std::vector<std::string> & fields = itr->second;
 
  141     std::map<std::string,Teuchos::RCP<const panzer::PureBasis> >::const_iterator found = bases.find(basisName);
 
  143                                "Could not find basis \""+basisName+
"\"!");
 
  147     std::vector<double> scalars(fields.size(),1.0); 
 
  148     for(std::size_t f=0;f<fields.size();f++) {
 
  149       std::unordered_map<std::string,double>::const_iterator f2s_itr = fieldToScalar_.find(fields[f]);
 
  153       if(f2s_itr!=fieldToScalar_.end()) {
 
  154         scalars[f] = f2s_itr->second;
 
  155         scaledFieldsHash.erase(fields[f]);
 
  163       std::string fields_concat = 
"";
 
  164       for(std::size_t f=0;f<fields.size();f++) {
 
  165         fields_concat += fields[f];
 
  170                                                       mesh_, basis, fields,scalars));
 
  173       this->
template registerEvaluator<EvalT>(fm, eval);
 
  174       fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
 
  183         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  187       std::string fields_concat = 
"";
 
  188       for(std::size_t f=0;f<fields.size();f++) {
 
  190         p.
set(
"Name",fields[f]);
 
  191         p.
set(
"Basis",basis);
 
  196         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  198         fields_concat += fields[f];
 
  205                                                                         mesh_,centroidRule,fields,scalars));
 
  207         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  208         fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]); 
 
  218         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  222       std::string fields_concat = 
"";
 
  223       for(std::size_t f=0;f<fields.size();f++) {
 
  225         p.
set(
"Name",fields[f]);
 
  226         p.
set(
"Basis",basis);
 
  231         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  233         fields_concat += fields[f];
 
  240                                                                         mesh_,centroidRule,fields,scalars));
 
  242         this->
template registerEvaluator<EvalT>(fm, evaluator);
 
  243         fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]); 
 
  252   for(std::unordered_set<std::string>::const_iterator itr=scaledFieldsHash.begin();
 
  253       itr!=scaledFieldsHash.end();itr++) { 
 
  254     out << 
"WARNING: STK Solution Writer did not scale the field \"" << *itr << 
"\" " 
  255         << 
"because it was not written." << std::endl;
 
  259 template <
typename EvalT>
 
  262                   std::map<std::string,std::vector<std::string> > & basisBucket)
 
  265    for(std::size_t i=0;i<providedDofs.size();i++) {
 
  266       std::string fieldName = providedDofs[i].first;
 
  269       basisBucket[basis->
name()].push_back(fieldName);
 
  273 template <
typename EvalT>
 
  277                          Kokkos::DynRankView<double,PHX::Device> & centroid)
 const 
  280    using Teuchos::rcp_dynamic_cast;
 
  282    centroid = Kokkos::DynRankView<double,PHX::Device>(
"centroid",1,baseDimension);
 
  286        itr!=bases.end();++itr) {
 
  291       Kokkos::DynRankView<double,PHX::Device> coords(
"coords",intrepidBasis->getCardinality(),
 
  292                  intrepidBasis->getBaseCellTopology().getDimension());
 
  293       intrepidBasis->getDofCoords(coords);
 
  297       for(
int i=0;i<coords.extent_int(0);i++)
 
  298          for(
int d=0;d<coords.extent_int(1);d++)
 
  299             centroid(0,d) += coords(i,d);
 
  302       for(
int d=0;d<coords.extent_int(1);d++)
 
  303          centroid(0,d) /= coords.extent(0);
 
  312 template <
typename EvalT>
 
  316   fieldToScalar_[fieldName] = fieldScalar;
 
  319 template <
typename EvalT>
 
  323   if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>())
 
  329 template <
typename EvalT>
 
  333   additionalFields_.push_back(std::make_pair(fieldName,basis));
 
  336 template <
typename EvalT>
 
  345   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 
 
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
virtual Teuchos::RCP< panzer::ResponseBase > buildResponseObject(const std::string &responseName) const 
 
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