Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl.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_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__
12 #define __Panzer_STK_ResponseEvaluatorFactory_SolutionWriter_impl_hpp__
13 
14 #include "Panzer_STK_Interface.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"
19 #include "Panzer_DOF.hpp"
21 
22 #include <unordered_set>
23 
24 namespace panzer_stk {
25 
26 namespace {
28  class Response_STKDummy : public panzer::ResponseBase {
29  public:
30  Response_STKDummy(const std::string & rn)
31  : ResponseBase(rn) {}
32  virtual void scatterResponse() {}
33  virtual void initializeResponse() {}
34  private:
35  Response_STKDummy();
36  Response_STKDummy(const Response_STKDummy &);
37  };
38 }
39 
40 template <typename EvalT>
42 buildResponseObject(const std::string & responseName) const
43 {
44  return Teuchos::rcp(new Response_STKDummy(responseName));
45 }
46 
47 template <typename EvalT>
49 buildAndRegisterEvaluators(const std::string& /* responseName */,
51  const panzer::PhysicsBlock & physicsBlock,
52  const Teuchos::ParameterList& /* user_data */) const
53 {
54  using Teuchos::RCP;
55  using Teuchos::rcp;
56 
57  typedef std::pair<std::string,RCP<const panzer::PureBasis> > StrConstPureBasisPair;
58 
59  // this will help so we can print out any unused scaled fields as a warning
60  std::unordered_set<std::string> scaledFieldsHash = scaledFieldsHash_;
61 
62  std::map<std::string,RCP<const panzer::PureBasis> > bases;
63  std::map<std::string,std::vector<std::string> > basisBucket;
64  {
65  const std::map<std::string,RCP<panzer::PureBasis> > & nc_bases = physicsBlock.getBases();
66  bases.insert(nc_bases.begin(),nc_bases.end());
67  }
68 
69  std::vector<StrConstPureBasisPair> allFields;
70 
71  // only add in solution fields if required
72 
73  if(!addCoordinateFields_ && addSolutionFields_) {
74  // inject all the fields, including the coordinates (we will remove them shortly)
75  allFields.insert(allFields.end(),physicsBlock.getProvidedDOFs().begin(),physicsBlock.getProvidedDOFs().end());
76 
77 
78  // get a list of strings with fields to remove
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]);
84 
85  // remove all coordinate fields
86  deleteRemovedFields(removedFields,allFields);
87  }
88  else if(addCoordinateFields_ && !addSolutionFields_) {
90  const std::vector<std::vector<std::string> > & coord_fields = physicsBlock.getCoordinateDOFs();
91 
92  // get the basis and field for each coordiante
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++) {
95  Teuchos::RCP<panzer::PureBasis> basis = // const_cast==yuck!
96  Teuchos::rcp_const_cast<panzer::PureBasis>(fieldLib->lookupBasis(coord_fields[c][d]));
97 
98  // make sure they are inserted in the allFields list
99  allFields.push_back(std::make_pair(coord_fields[c][d],basis));
100  }
101  }
102  }
103  else if(addSolutionFields_)
104  allFields.insert(allFields.end(),physicsBlock.getProvidedDOFs().begin(),physicsBlock.getProvidedDOFs().end());
105 
106  // Add in tangent fields
107  if(addSolutionFields_)
108  allFields.insert(allFields.end(),physicsBlock.getTangentFields().begin(),physicsBlock.getTangentFields().end());
109 
110  // add in bases for any addtional fields
111  for(std::size_t i=0;i<additionalFields_.size();i++)
112  bases[additionalFields_[i].second->name()] = additionalFields_[i].second;
113 
114  allFields.insert(allFields.end(),additionalFields_.begin(),additionalFields_.end());
115 
116  deleteRemovedFields(removedFields_,allFields);
117 
118  bucketByBasisType(allFields,basisBucket);
119 
120  // add this for HCURL and HDIV basis, only want to add them once: evaluate vector fields at centroid
122  RCP<panzer::PointRule> centroidRule;
123  for(std::map<std::string,Teuchos::RCP<const panzer::PureBasis> >::const_iterator itr=bases.begin();
124  itr!=bases.end();++itr) {
125 
126  if(itr->second->isVectorBasis()) {
127  centroidRule = rcp(new panzer::PointRule("Centroid",1,physicsBlock.cellData()));
128 
129  // compute centroid
130  Kokkos::DynRankView<double,PHX::Device> centroid;
131  computeReferenceCentroid(bases,physicsBlock.cellData().baseCellDimension(),centroid);
132 
133  // build pointe values evaluator
135  rcp(new panzer::PointValues_Evaluator<EvalT,panzer::Traits>(centroidRule,centroid));
136  this->template registerEvaluator<EvalT>(fm, evaluator);
137 
138  break; // get out of the loop, only need one evaluator
139  }
140  }
141 
142  // add evaluators for each field
144 
145  for(std::map<std::string,std::vector<std::string> >::const_iterator itr=basisBucket.begin();
146  itr!=basisBucket.end();++itr) {
147 
148  std::string basisName = itr->first;
149  const std::vector<std::string> & fields = itr->second;
150 
151  std::map<std::string,Teuchos::RCP<const panzer::PureBasis> >::const_iterator found = bases.find(basisName);
152  TEUCHOS_TEST_FOR_EXCEPTION(found==bases.end(),std::logic_error,
153  "Could not find basis \""+basisName+"\"!");
154  Teuchos::RCP<const panzer::PureBasis> basis = found->second;
155 
156  // determine if user has modified field scalar for each field to be written to STK
157  std::vector<double> scalars(fields.size(),1.0); // fill with 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]);
160 
161  // if scalar is found, include it in the vector and remove the field from the
162  // hash table so it won't be included in the warning message.
163  if(f2s_itr!=fieldToScalar_.end()) {
164  scalars[f] = f2s_itr->second;
165  scaledFieldsHash.erase(fields[f]);
166  }
167  }
168 
169  // write out nodal fields
170  if(basis->getElementSpace()==panzer::PureBasis::HGRAD ||
171  basis->getElementSpace()==panzer::PureBasis::CONST) {
172 
173  std::string fields_concat = "";
174  for(std::size_t f=0;f<fields.size();f++) {
175  fields_concat += fields[f];
176  }
177 
179  Teuchos::rcp(new ScatterFields<EvalT,panzer::Traits>("STK HGRAD Scatter Basis " +basis->name()+": "+fields_concat,
180  mesh_, basis, fields,scalars));
181 
182  // register and require evaluator fields
183  this->template registerEvaluator<EvalT>(fm, eval);
184  fm.template requireField<EvalT>(*eval->evaluatedFields()[0]);
185  }
186  else if(basis->getElementSpace()==panzer::PureBasis::HCURL) {
187  TEUCHOS_ASSERT(centroidRule!=Teuchos::null);
188 
189  // register basis values evaluator
190  {
193  this->template registerEvaluator<EvalT>(fm, evaluator);
194  }
195 
196  // add a DOF_PointValues for each field
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);
202  p.set("Point Rule",centroidRule.getConst());
205 
206  this->template registerEvaluator<EvalT>(fm, evaluator);
207 
208  fields_concat += fields[f];
209  }
210 
211  // add the scatter field evaluator for this basis
212  {
214  = Teuchos::rcp(new ScatterVectorFields<EvalT,panzer::Traits>("STK HCURL Scatter Basis " +basis->name()+": "+fields_concat,
215  mesh_,centroidRule,fields,scalars));
216 
217  this->template registerEvaluator<EvalT>(fm, evaluator);
218  fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]); // require the dummy evaluator
219  }
220  }
221  else if(basis->getElementSpace()==panzer::PureBasis::HDIV) {
222  TEUCHOS_ASSERT(centroidRule!=Teuchos::null);
223 
224  // register basis values evaluator
225  {
228  this->template registerEvaluator<EvalT>(fm, evaluator);
229  }
230 
231  // add a DOF_PointValues for each field
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);
237  p.set("Point Rule",centroidRule.getConst());
240 
241  this->template registerEvaluator<EvalT>(fm, evaluator);
242 
243  fields_concat += fields[f];
244  }
245 
246  // add the scatter field evaluator for this basis
247  {
249  = Teuchos::rcp(new ScatterVectorFields<EvalT,panzer::Traits>("STK HDIV Scatter Basis " +basis->name()+": "+fields_concat,
250  mesh_,centroidRule,fields,scalars));
251 
252  this->template registerEvaluator<EvalT>(fm, evaluator);
253  fm.template requireField<EvalT>(*evaluator->evaluatedFields()[0]); // require the dummy evaluator
254  }
255  }
256  }
257 
258  // print warning message for any unused scaled fields
259  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
260  out.setOutputToRootOnly(0);
261 
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;
266  }
267 }
268 
269 template <typename EvalT>
271 bucketByBasisType(const std::vector<std::pair<std::string,Teuchos::RCP<const panzer::PureBasis> > > & providedDofs,
272  std::map<std::string,std::vector<std::string> > & basisBucket)
273 {
274  // this should be self explanatory
275  for(std::size_t i=0;i<providedDofs.size();i++) {
276  std::string fieldName = providedDofs[i].first;
277  Teuchos::RCP<const panzer::PureBasis> basis = providedDofs[i].second;
278 
279  basisBucket[basis->name()].push_back(fieldName);
280  }
281 }
282 
283 template <typename EvalT>
286  int baseDimension,
287  Kokkos::DynRankView<double,PHX::Device> & centroid) const
288 {
289  using Teuchos::RCP;
290  using Teuchos::rcp_dynamic_cast;
291 
292  centroid = Kokkos::DynRankView<double,PHX::Device>("centroid",1,baseDimension);
293  auto l_centroid = centroid;
294 
295  // loop over each possible basis
296  for(std::map<std::string,RCP<const panzer::PureBasis> >::const_iterator itr=bases.begin();
297  itr!=bases.end();++itr) {
298 
299  RCP<Intrepid2::Basis<PHX::exec_space,double,double>> intrepidBasis = itr->second->getIntrepid2Basis();
300 
301  // we've got coordinates, lets commpute the "centroid"
302  Kokkos::DynRankView<double,PHX::Device> coords("coords",intrepidBasis->getCardinality(),
303  intrepidBasis->getBaseCellTopology().getDimension());
304  intrepidBasis->getDofCoords(coords);
305  TEUCHOS_ASSERT(coords.rank()==2);
306  TEUCHOS_ASSERT(coords.extent_int(1)==baseDimension);
307 
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);
311  });
312  Kokkos::fence();
313 
314  return;
315  }
316 
317  // no centroid was found...die
318  TEUCHOS_ASSERT(false);
319 }
320 
321 template <typename EvalT>
323 scaleField(const std::string & fieldName,double fieldScalar)
324 {
325  fieldToScalar_[fieldName] = fieldScalar;
326  scaledFieldsHash_.insert(fieldName);
327 }
328 
329 template <typename EvalT>
332 {
333  if(PHX::print<EvalT>()==PHX::print<panzer::Traits::Residual>())
334  return true;
335 
336  return false;
337 }
338 
339 template <typename EvalT>
341 addAdditionalField(const std::string & fieldName,const Teuchos::RCP<const panzer::PureBasis> & basis)
342 {
343  additionalFields_.push_back(std::make_pair(fieldName,basis));
344 }
345 
346 template <typename EvalT>
348 deleteRemovedFields(const std::vector<std::string> & removedFields,
349  std::vector<std::pair<std::string,Teuchos::RCP<const panzer::PureBasis> > > & fields) const
350 {
352  functor.removedFields_ = removedFields;
353 
354  // This is the Erase-Remove Idiom: see http://en.wikipedia.org/wiki/Erase-remove_idiom
355  fields.erase(std::remove_if(fields.begin(),fields.end(),functor),fields.end());
356 }
357 
358 }
359 
360 #endif
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.
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
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.
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 > &centroid) const