Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_EquationSet_DefaultImpl_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_EQUATIONSET_DEFAULT_IMPL_IMPL_HPP
44 #define PANZER_EQUATIONSET_DEFAULT_IMPL_IMPL_HPP
45 
46 #include "Panzer_DOF.hpp"
48 #include "Panzer_DOFGradient.hpp"
49 #include "Panzer_DOFCurl.hpp"
50 #include "Panzer_DOFDiv.hpp"
51 #include "Panzer_GatherBasisCoordinates.hpp"
52 #include "Panzer_GatherIntegrationCoordinates.hpp"
53 #include "Panzer_GatherOrientation.hpp"
54 #include "Panzer_GlobalIndexer.hpp"
55 
56 #include "Phalanx_MDField.hpp"
57 #include "Phalanx_DataLayout.hpp"
58 #include "Phalanx_DataLayout_MDALayout.hpp"
60 
61 // For convenience we automate some evalautor registration
62 #include "Panzer_Sum.hpp"
63 
64 // ***********************************************************************
65 template <typename EvalT>
68  const int& default_integration_order,
69  const panzer::CellData& cell_data,
70  const Teuchos::RCP<panzer::GlobalData>& global_data,
71  const bool build_transient_support) :
72  panzer::GlobalDataAcceptorDefaultImpl(global_data),
73  m_input_params(params),
74  m_default_integration_order(default_integration_order),
75  m_cell_data(cell_data),
76  m_build_transient_support(build_transient_support)
77 {
79 
81 }
82 
83 // ***********************************************************************
84 template <typename EvalT>
87 {
88  // Get defaultparameters
89  m_type = m_input_params->get<std::string>("Type");
90 
91  // for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
92  // itr!=m_provided_dofs_desc.end();++itr) {
93  // itr->second.print(std::cout); std::cout << std::endl;
94  // }
95 
96  this->m_provided_dofs.clear();
97  this->m_int_rules.clear();
98 
99  // load up the provided dofs and unique int rules from the descriptor map
100  for(typename std::map<std::string,DOFDescriptor>::iterator itr=m_provided_dofs_desc.begin();
101  itr!=m_provided_dofs_desc.end();++itr) {
102 
103  // Create the bases
104  TEUCHOS_ASSERT(nonnull(itr->second.basis));
105  this->m_provided_dofs.push_back(std::make_pair(itr->first, itr->second.basis));
106 
107  //{
108  // Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
109  // out.setOutputToRootOnly(0);
110  // itr->second.print(out);
111  // out << "Int Order = " << itr->second.integrationOrder << " (" << itr->second.intRule->order() << ")" << std::endl;
112  //}
113 
114  // Create the unique integration rule map and complete descriptor objects
115  TEUCHOS_ASSERT(nonnull(itr->second.intRule));
116  m_int_rules[itr->second.intRule->order()] = itr->second.intRule;
117 
118  }
119 
120  // Setup the basis to dof mapping
121  for (DescriptorIterator dof_iter = m_provided_dofs_desc.begin(); dof_iter != m_provided_dofs_desc.end(); ++dof_iter) {
122 
123  std::string basis_name = dof_iter->second.basis->name();
124  Teuchos::RCP<panzer::PureBasis> basis = dof_iter->second.basis;
125  std::string dof_name = dof_iter->first;
126 
127  if (is_null(m_basis_to_dofs[basis_name].first)) {
128  m_basis_to_dofs[basis_name].first = basis;
129  m_basis_to_dofs[basis_name].second = Teuchos::rcp(new std::vector<std::string>);
130  }
131 
132  m_basis_to_dofs[basis_name].second->push_back(dof_name);
133  }
134 
135  // Generate a unique list of bases
136  for (DescriptorIterator dof_iter = m_provided_dofs_desc.begin(); dof_iter != m_provided_dofs_desc.end(); ++dof_iter) {
137  m_unique_bases[dof_iter->second.basis->name()] = dof_iter->second.basis;
138  }
139 
140  // Setup the default parameter list for closure models
141  this->m_eval_plist->set("Block ID", getElementBlockId());
142  this->setupDeprecatedDOFsSupport();
143 }
144 
145 // ***********************************************************************
146 template <typename EvalT>
148 {
149  TEUCHOS_ASSERT(m_provided_dofs.size() > 0);
150  TEUCHOS_ASSERT(m_int_rules.size() > 0);
151 
152  // Deprecated support assumes all equations in set use the same
153  // basis and integration rule
154  Teuchos::RCP<panzer::PureBasis> pure_basis = m_provided_dofs.begin()->second;
155  Teuchos::RCP<panzer::IntegrationRule> int_rule = m_int_rules.begin()->second;
156  Teuchos::RCP<panzer::BasisIRLayout> basis = panzer::basisIRLayout(pure_basis,*int_rule);
157 
158  this->m_eval_plist->set("Basis", basis);
159  this->m_eval_plist->set("IR", int_rule);
160 }
161 
162 // ***********************************************************************
163 template <typename EvalT>
166  const panzer::FieldLibrary& /* fl */,
168  const Teuchos::ParameterList& /* user_data */) const
169 {
171  using Teuchos::RCP;
172  using Teuchos::rcp;
173 
174  // ********************
175  // DOFs (unknowns)
176  // ********************
177 
178  // Gather, includes construction of orientation gathers
179  for (BasisIterator basis_it = m_basis_to_dofs.begin(); basis_it != m_basis_to_dofs.end(); ++basis_it) {
180 
181  // Set tangent field names (first dimension is DOF, second is parameter)
183  if (m_tangent_param_names.size() > 0) {
184  tangent_field_names = rcp(new std::vector< std::vector<std::string> >(basis_it->second.second->size()));
185  for (std::size_t i=0; i<basis_it->second.second->size(); ++i) {
186  for (std::size_t j=0; j<m_tangent_param_names.size(); ++j) {
187  const std::string tname =
188  (*(basis_it->second.second))[i] + " SENSITIVITY " + m_tangent_param_names[j];
189  (*tangent_field_names)[i].push_back(tname);
190  }
191  }
192  }
193 
194  {
195  ParameterList p("Gather");
196  p.set("Basis", basis_it->second.first);
197  p.set("DOF Names", basis_it->second.second);
198  p.set("Indexer Names", basis_it->second.second);
199  p.set("Sensitivities Name", "");
200  p.set("First Sensitivities Available", true);
201  p.set("Second Sensitivities Available", true);
202 
203  // Set tangent field names
204  if (tangent_field_names != Teuchos::null)
205  p.set("Tangent Names", tangent_field_names);
206 
208 
209  this->template registerEvaluator<EvalT>(fm, op);
210  }
211 
212  // Create a second gather evaluator for each tangent field,
213  // we never compute derivatives with respect to this field
214  if (tangent_field_names != Teuchos::null) {
215  for (std::size_t i=0; i<m_tangent_param_names.size(); ++i) {
216 
217  Teuchos::RCP< std::vector<std::string> > names = rcp(new std::vector<std::string>);
218  for (std::size_t j=0; j<basis_it->second.second->size(); ++j)
219  names->push_back((*tangent_field_names)[j][i]);
220 
221  ParameterList p(std::string("Gather Tangent ") + this->m_tangent_param_names[i]);
222  p.set("Basis", basis_it->second.first);
223  p.set("DOF Names", names);
224  p.set("Indexer Names", basis_it->second.second);
225  p.set("Sensitivities Name", "");
226  p.set("First Sensitivities Available", false);
227  p.set("Second Sensitivities Available", false);
228  p.set("Global Data Key", "X TANGENT GATHER CONTAINER: " + this->m_tangent_param_names[i]);
229 
231 
232  this->template registerEvaluator<EvalT>(fm, op);
233  }
234  }
235  }
236 
237  // **************************
238  // Coordinates for integration points and basis functions
239  // **************************
240  {
241  // add basis coordinates
242  for (std::map<std::string,Teuchos::RCP<panzer::PureBasis> >::const_iterator basis = m_unique_bases.begin();
243  basis != m_unique_bases.end(); ++ basis) {
246  this->template registerEvaluator<EvalT>(fm, basis_op);
247  }
248 
249  // add integration coordinates
250  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir = m_int_rules.begin();
251  ir != m_int_rules.end(); ++ir) {
254  this->template registerEvaluator<EvalT>(fm, quad_op);
255  }
256 
257  // NOTE: You can look up the name of either coordinate field name by doing
258  // GatherBasisCoordinates<EvalT,Traits>::fieldName();
259  // GatherIntegrationCoordinates<EvalT,Traits>::fieldName();
260  }
261 
262  // **************************
263  // Time derivative terms
264  // **************************
265 
266  // Gather of time derivative terms: One evaluator for each unique basis
267  for (BasisIterator basis_it = m_basis_to_dofs.begin(); basis_it != m_basis_to_dofs.end(); ++basis_it) {
268 
269  RCP< std::vector<std::string> > t_dof_names = rcp(new std::vector<std::string>); // time derivative indexer names
270  RCP< std::vector<std::string> > t_field_names = rcp(new std::vector<std::string>); // time derivative field names
271  RCP< std::vector< std::vector<std::string> > > tangent_field_names = rcp(new std::vector< std::vector<std::string> >); // tangent field names
272 
273  // determine which fields associated with this basis need time derivatives
274  for (typename std::vector<std::string>::const_iterator dof_name = basis_it->second.second->begin();
275  dof_name != basis_it->second.second->end(); ++dof_name) {
276 
277  DescriptorIterator desc = m_provided_dofs_desc.find(*dof_name);
278  TEUCHOS_ASSERT(desc != m_provided_dofs_desc.end());
279 
280  // does this field need a time derivative?
281  if(desc->second.timeDerivative.first) {
282  // time derivative needed
283  t_dof_names->push_back(*dof_name);
284  t_field_names->push_back(desc->second.timeDerivative.second);
285 
286  // Set tangent field names (first dimension is DOF, second is parameter)
287  if (m_tangent_param_names.size() > 0) {
288  std::vector<std::string> tfn;
289  for (std::size_t j=0; j<m_tangent_param_names.size(); ++j) {
290  const std::string tname =
291  desc->second.timeDerivative.second + " SENSITIVITY " + m_tangent_param_names[j];
292  tfn.push_back(tname);
293  }
294  tangent_field_names->push_back(tfn);
295  }
296  }
297  }
298 
299  {
300  ParameterList p("Gather");
301  p.set("Basis", basis_it->second.first);
302  p.set("DOF Names", t_field_names);
303  p.set("Indexer Names", t_dof_names);
304  p.set("Use Time Derivative Solution Vector", true);
305 
306  // Set tangent field names
307  if (m_tangent_param_names.size() > 0)
308  p.set("Tangent Names", tangent_field_names);
309 
311 
312  this->template registerEvaluator<EvalT>(fm, op);
313  }
314 
315  // Create a second gather evaluator for each tangent field,
316  // we never compute derivatives with respect to this field
317  if (m_tangent_param_names.size() > 0) {
318  for (std::size_t i=0; i<m_tangent_param_names.size(); ++i) {
319 
321  rcp(new std::vector<std::string>);
322  for (std::size_t j=0; j<tangent_field_names->size(); ++j)
323  names->push_back((*tangent_field_names)[j][i]);
324 
325  ParameterList p(std::string("Gather Tangent ") + this->m_tangent_param_names[i]);
326  p.set("Basis", basis_it->second.first);
327  p.set("DOF Names", names);
328  p.set("Indexer Names", t_dof_names);
329  p.set("Use Time Derivative Solution Vector", true);
330  p.set("Global Data Key", "DXDT TANGENT GATHER CONTAINER: " + this->m_tangent_param_names[i]);
331 
333 
334  this->template registerEvaluator<EvalT>(fm, op);
335  }
336  }
337  }
338 
339  // **************************
340  // Orientation terms
341  // **************************
342 
343  for (BasisIterator basis_it = m_basis_to_dofs.begin(); basis_it != m_basis_to_dofs.end(); ++basis_it) {
344  if(basis_it->second.first->requiresOrientations()) {
345  ParameterList p("Gather Orientation");
346  p.set("Basis", basis_it->second.first);
347  p.set("DOF Names", basis_it->second.second);
348  p.set("Indexer Names", basis_it->second.second);
349 
351 
352  this->template registerEvaluator<EvalT>(fm, op);
353  }
354  }
355 
356 }
357 
358 // ***********************************************************************
359 template <typename EvalT>
362  const panzer::FieldLayoutLibrary& fl,
365  const Teuchos::ParameterList& /* user_data */) const
366 {
368  using Teuchos::RCP;
369  using Teuchos::rcp;
370 
372  if(lof!=Teuchos::null)
373  globalIndexer = lof->getRangeGlobalIndexer();
374 
375  // DOFs: Scalar value @ basis --> Scalar value @ IP
376  for (DescriptorIterator dof_iter = m_provided_dofs_desc.begin(); dof_iter != m_provided_dofs_desc.end(); ++dof_iter) {
377 
378  ParameterList p;
379  p.set("Name", dof_iter->first);
380  p.set("Basis", fl.lookupLayout(dof_iter->first));
381  p.set("IR", ir);
382 
383  if(globalIndexer!=Teuchos::null) {
384  // build the offsets for this field
385  int fieldNum = globalIndexer->getFieldNum(dof_iter->first);
387  rcp(new std::vector<int>(globalIndexer->getGIDFieldOffsets(m_block_id,fieldNum)));
388  p.set("Jacobian Offsets Vector", offsets);
389  }
390  // else default to the slow DOF call
391 
394 
395  this->template registerEvaluator<EvalT>(fm, op);
396  }
397 
398  // Gradients of DOFs: Scalar value @ basis --> Vector value @ IP
399 
400  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
401  itr!=m_provided_dofs_desc.end();++itr) {
402 
403  if(itr->second.basis->supportsGrad()) {
404 
405  // is gradient required for this variable
406  if(!itr->second.grad.first)
407  continue; // its not required, quit the loop
408 
409  const std::string dof_name = itr->first;
410  const std::string dof_grad_name = itr->second.grad.second;
411 
412  ParameterList p;
413  p.set("Name", dof_name);
414  p.set("Gradient Name", dof_grad_name);
415  p.set("Basis", fl.lookupLayout(dof_name));
416  p.set("IR", ir);
417 
420 
421  this->template registerEvaluator<EvalT>(fm, op);
422  }
423  }
424 
425  // Curl of DOFs: Vector value @ basis --> Vector value @ IP (3D) or Scalar value @ IP (2D)
426 
427  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
428  itr!=m_provided_dofs_desc.end();++itr) {
429 
430  if(itr->second.basis->supportsCurl()) {
431 
432  // is curl required for this variable
433  if(!itr->second.curl.first)
434  continue; // its not required, quit the loop
435 
436  const std::string dof_name = itr->first;
437  const std::string dof_curl_name = itr->second.curl.second;
438 
439  ParameterList p;
440  p.set("Name", dof_name);
441  p.set("Curl Name", dof_curl_name);
442  p.set("Basis", fl.lookupLayout(dof_name));
443  p.set("IR", ir);
444 
445  // this will help accelerate the DOFCurl evaluator when Jacobians are needed
446  if(globalIndexer!=Teuchos::null) {
447  // build the offsets for this field
448  int fieldNum = globalIndexer->getFieldNum(dof_name);
450  rcp(new std::vector<int>(globalIndexer->getGIDFieldOffsets(m_block_id,fieldNum)));
451  p.set("Jacobian Offsets Vector", offsets);
452  }
453  // else default to the slow DOF call
454 
455 
458 
459  this->template registerEvaluator<EvalT>(fm, op);
460  }
461 
462  }
463 
464  // Div of DOFs: Vector value @ basis --> Scalar value @ IP
465 
466  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
467  itr!=m_provided_dofs_desc.end();++itr) {
468 
469  if(itr->second.basis->supportsDiv()) {
470 
471  // is div required for this variable
472  if(!itr->second.div.first)
473  continue; // its not required, quit the loop
474 
475  const std::string dof_name = itr->first;
476  const std::string dof_div_name = itr->second.div.second;
477 
478  ParameterList p;
479  p.set("Name", dof_name);
480  p.set("Div Name", dof_div_name);
481  p.set("Basis", fl.lookupLayout(dof_name));
482  p.set("IR", ir);
483 
484  // this will help accelerate the DOFDiv evaluator when Jacobians are needed
485  if(globalIndexer!=Teuchos::null) {
486  // build the offsets for this field
487  int fieldNum = globalIndexer->getFieldNum(dof_name);
489  rcp(new std::vector<int>(globalIndexer->getGIDFieldOffsets(m_block_id,fieldNum)));
490  p.set("Jacobian Offsets Vector", offsets);
491  }
492  // else default to the slow DOF call
493 
494 
497 
498  this->template registerEvaluator<EvalT>(fm, op);
499  }
500  }
501 
502  // Time derivative of DOFs: Scalar value @ basis --> Scalar value @ IP
503  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
504  itr!=m_provided_dofs_desc.end();++itr) {
505  // is td required for this variable
506  if(!itr->second.timeDerivative.first)
507  continue; // its not required, quit the loop
508 
509  const std::string td_name = itr->second.timeDerivative.second;
510 
511  ParameterList p;
512  p.set("Name", td_name);
513  p.set("Basis", fl.lookupLayout(itr->first));
514  p.set("IR", ir);
515 
516  if(globalIndexer!=Teuchos::null) {
517  // build the offsets for this field
518  int fieldNum = globalIndexer->getFieldNum(itr->first);
520  rcp(new std::vector<int>(globalIndexer->getGIDFieldOffsets(m_block_id,fieldNum)));
521  p.set("Jacobian Offsets Vector", offsets);
522  }
523  // else default to the slow DOF call
524 
525  // set the orientiation field name explicitly if orientations are
526  // required for the basis
527  if(itr->second.basis->requiresOrientations())
528  p.set("Orientation Field Name", itr->first+" Orientation");
529 
532 
533  this->template registerEvaluator<EvalT>(fm, op);
534  }
535 
536 }
537 
538 // ***********************************************************************
539 template <typename EvalT>
542  const panzer::FieldLibrary& /* fl */,
544  const Teuchos::ParameterList& user_data) const
545 {
547  using Teuchos::RCP;
548  using Teuchos::rcp;
549 
550  // this turns off the scatter contribution, and does
551  // only the gather
552  bool ignoreScatter = false;
553  if(user_data.isParameter("Ignore Scatter"))
554  ignoreScatter = user_data.get<bool>("Ignore Scatter");
555 
556  if(!ignoreScatter) {
557 
558  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
559  itr!=m_provided_dofs_desc.end();++itr) {
560 
561  RCP<std::map<std::string,std::string> > names_map = rcp(new std::map<std::string,std::string>);
562  RCP< std::vector<std::string> > residual_names = rcp(new std::vector<std::string>);
563 
564  // sanity check to make sure a residual name was registered for each provided variable
565  TEUCHOS_ASSERT(itr->second.residualName.first);
566 
567  names_map->insert(std::make_pair(itr->second.residualName.second,itr->first));
568  residual_names->push_back(itr->second.residualName.second);
569 
570  {
571  ParameterList p("Scatter");
572  p.set("Scatter Name", itr->second.scatterName);
573  p.set("Basis", itr->second.basis.getConst());
574  p.set("Dependent Names", residual_names);
575  p.set("Dependent Map", names_map);
576 
578 
579  this->template registerEvaluator<EvalT>(fm, op);
580  }
581 
582  // Require variables
583  {
584  PHX::Tag<typename EvalT::ScalarT> tag(itr->second.scatterName,
586  fm.template requireField<EvalT>(tag);
587  }
588 
589  }
590 
591  }
592 
593 }
594 
595 // ***********************************************************************
596 
597 template <typename EvalT>
600  const panzer::FieldLayoutLibrary& fl,
603  const Teuchos::ParameterList& models,
604  const Teuchos::ParameterList& user_data) const
605 {
606  for (std::vector<std::string>::const_iterator model_name = m_closure_model_ids.begin();
607  model_name != m_closure_model_ids.end(); ++model_name) {
608 
609  this->buildAndRegisterClosureModelEvaluators(fm,fl,ir,factory,*model_name,models,user_data);
610  }
611 }
612 
613 // ***********************************************************************
614 
615 template <typename EvalT>
618  const panzer::FieldLayoutLibrary& fl,
621  const std::string& model_name,
622  const Teuchos::ParameterList& models,
623  const Teuchos::ParameterList& user_data) const
624 {
626  factory.getAsObject<EvalT>()->buildClosureModels(model_name,
627  models,
628  fl,
629  ir,
630  *(this->m_eval_plist),
631  user_data,
632  this->getGlobalData(),
633  fm);
634 
635  for (std::vector< Teuchos::RCP<PHX::Evaluator<panzer::Traits> > >::size_type i=0; i < evaluators->size(); ++i)
636  this->template registerEvaluator<EvalT>(fm, (*evaluators)[i]);
637 }
638 
639 // ***********************************************************************
640 template <typename EvalT>
643  const panzer::FieldLibrary& fl,
645  const std::string& model_name,
646  const Teuchos::ParameterList& models,
648  const Teuchos::ParameterList& user_data) const
649 {
650  // add basis coordinates
651  for (std::map<std::string,Teuchos::RCP<panzer::PureBasis> >::const_iterator basis = m_unique_bases.begin();
652  basis != m_unique_bases.end(); ++ basis) {
655  this->template registerEvaluator<EvalT>(fm, basis_op);
656  }
657 
658  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
659  itr!=m_provided_dofs_desc.end();++itr) {
660 
661  Teuchos::ParameterList p("Scatter");
662  p.set("Scatter Name", itr->second.scatterName);
663  p.set("Basis", itr->second.basis.getConst());
664  Teuchos::RCP<std::vector<std::string> > name = Teuchos::rcp(new std::vector<std::string>);
665  name->push_back(itr->first);
666  p.set("Dependent Names", name);
667 
668  // Create an identity map
669  Teuchos::RCP<std::map<std::string,std::string> > names_map = Teuchos::rcp(new std::map<std::string,std::string>);
670  names_map->insert(std::make_pair(itr->first,itr->first));
671  p.set("Dependent Map", names_map);
672 
673  // Set flag for ScatterDirichlet evaluators
674  p.set("Scatter Initial Condition", true);
675 
676  // Use ScatterDirichlet to scatter the initial condition
678 
679  this->template registerEvaluator<EvalT>(fm, op);
680 
681 
682  // Require field
683  PHX::Tag<typename EvalT::ScalarT> tag(itr->second.scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0)));
684  fm.template requireField<EvalT>(tag);
685  }
686 
687  // Add in closure models. This is a hack that we should revisit.
688  {
689  // Closure models are normally evaluated at integration points,
690  // but some evaluator models support evaluation at both basis and
691  // integration points. For initial guesses, we should only
692  // evaluate at basis points, so integration rule is meaningless.
693  // We use this to build all closure model evaluators in model
694  // (including integration point based ones that will never be
695  // used). In the future we may need ir for using L2 projection to
696  // basis points for initial guesses (for non-nodal bases).
698  if (m_int_rules.size() > 0)
699  dummy_ir = m_int_rules.begin()->second;
701 
703  factory.getAsObject<EvalT>()->buildClosureModels(model_name, models, *fll, dummy_ir, *(this->m_eval_plist), user_data, this->getGlobalData(), fm);
704 
705  for (std::vector< Teuchos::RCP<PHX::Evaluator<panzer::Traits> > >::size_type i=0; i < evaluators->size(); ++i)
706  this->template registerEvaluator<EvalT>(fm, (*evaluators)[i]);
707  }
708 
709  // **************************
710  // Add Orientation terms
711  // **************************
712 
713  for (BasisIterator basis_it = m_basis_to_dofs.begin(); basis_it != m_basis_to_dofs.end(); ++basis_it) {
714  if(basis_it->second.first->requiresOrientations()) {
715  Teuchos::ParameterList p("Gather Orientation");
716  p.set("Basis", basis_it->second.first);
717  p.set("DOF Names", basis_it->second.second);
718  p.set("Indexer Names", basis_it->second.second);
719 
721 
722  this->template registerEvaluator<EvalT>(fm, op);
723  }
724  }
725 
726 }
727 
728 // ***********************************************************************
729 template <typename EvalT>
732 {
733  return m_eval_plist;
734 }
735 
736 // ***********************************************************************
737 template <typename EvalT>
738 const std::vector<std::pair<std::string,Teuchos::RCP<panzer::PureBasis> > >&
740 {
741  return m_provided_dofs;
742 }
743 
744 // ***********************************************************************
745 template <typename EvalT>
746 const std::vector<std::vector<std::string> > &
748 {
749  return m_coordinate_dofs;
750 }
751 
752 // ***********************************************************************
753 template <typename EvalT>
754 const std::map<int,Teuchos::RCP<panzer::IntegrationRule> > &
756 {
757  return m_int_rules;
758 }
759 
760 // ***********************************************************************
761 template <typename EvalT>
763 setElementBlockId(const std::string & blockId)
764 {
765  TEUCHOS_ASSERT(m_block_id=="");
766  m_block_id = blockId;
767  this->m_eval_plist->set("Block ID", getElementBlockId()); // set the value in parameter list
768  // used by closure model factory
769 }
770 
771 // ***********************************************************************
772 template <typename EvalT>
775 {
776  return m_block_id;
777 }
778 
779 // ***********************************************************************
780 template <typename EvalT>
782 {
783  return m_type;
784 }
785 
786 // ***********************************************************************
787 template <typename EvalT>
788 void panzer::EquationSet_DefaultImpl<EvalT>::setTangentParamNames(const std::vector<std::string>& tangent_param_names)
789 {
790  m_tangent_param_names = tangent_param_names;
791 }
792 
793 // ***********************************************************************
794 template <typename EvalT>
796 {
797  return m_build_transient_support;
798 }
799 
800 // ***********************************************************************
801 template <typename EvalT>
803 getAddedDOFs(std::vector<std::string> & dofNames) const
804 {
805  dofNames.clear();
806  for(typename std::map<std::string,DOFDescriptor>::const_iterator itr=m_provided_dofs_desc.begin();
807  itr!=m_provided_dofs_desc.end();++itr)
808  dofNames.push_back(itr->first);
809 }
810 
811 // ***********************************************************************
812 template <typename EvalT>
814 updateDOF(const std::string & dofName,
815  int basisOrder,
816  int integrationOrder)
817 {
818  typename std::map<std::string,DOFDescriptor>::const_iterator itr = m_provided_dofs_desc.find(dofName);
819 
820  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
821  "EquationSet_DefaultImpl::updateDOF: DOF \"" << dofName << "\" has not been specified "
822  "by derived equation set \"" << this->getType() << "\".");
823 
824  // allocate and populate a dof descriptor associated with the field "dofName"
825  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
826  desc.basisOrder = basisOrder;
827  desc.basis = Teuchos::rcp(new panzer::PureBasis(desc.basisType,desc.basisOrder,m_cell_data));
828 
829  if (integrationOrder == -1)
830  desc.integrationOrder = m_default_integration_order;
831  else
832  desc.integrationOrder = integrationOrder;
833 
834  desc.intRule = Teuchos::rcp(new panzer::IntegrationRule(desc.integrationOrder,m_cell_data));
835 }
836 
837 // ***********************************************************************
838 template <typename EvalT>
840 getBasisOrder(const std::string & dofName) const
841 {
842  typename std::map<std::string,DOFDescriptor>::const_iterator itr = m_provided_dofs_desc.find(dofName);
843 
844  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
845  "EquationSet_DefaultImpl::getBasisOrder: DOF \"" << dofName << "\" has not been specified "
846  "by derived equation set \"" << this->getType() << "\".");
847 
848  return itr->second.basisOrder;
849 }
850 
851 // ***********************************************************************
852 template <typename EvalT>
854 getIntegrationOrder(const std::string & dofName) const
855 {
856  typename std::map<std::string,DOFDescriptor>::const_iterator itr = m_provided_dofs_desc.find(dofName);
857 
858  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
859  "EquationSet_DefaultImpl::getIntegrationOrder: DOF \"" << dofName << "\" has not been specified "
860  "by derived equation set \"" << this->getType() << "\".");
861 
862  return itr->second.integrationOrder;
863 }
864 
865 // ***********************************************************************
866 template <typename EvalT>
868 addDOF(const std::string & dofName,
869  const std::string & basisType,
870  const int & basisOrder,
871  const int integrationOrder,
872  const std::string residualName,
873  const std::string scatterName)
874 {
875  typename std::map<std::string,DOFDescriptor>::const_iterator itr = m_provided_dofs_desc.find(dofName);
876 
877  TEUCHOS_TEST_FOR_EXCEPTION(itr!=m_provided_dofs_desc.end(),std::runtime_error,
878  "EquationSet_DefaultImpl::addProvidedDOF: DOF \"" << dofName << "\" was previously specified "
879  "by derived equation set \"" << this->getType() << "\".");
880 
881  // allocate and populate a dof descriptor associated with the field "dofName"
882  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
883  desc.dofName = dofName;
884  desc.basisType = basisType;
885  desc.basisOrder = basisOrder;
886  desc.basis = Teuchos::rcp(new panzer::PureBasis(desc.basisType,desc.basisOrder,m_cell_data));
887 
888  if (integrationOrder == -1)
889  desc.integrationOrder = m_default_integration_order;
890  else
891  desc.integrationOrder = integrationOrder;
892 
893  desc.intRule = Teuchos::rcp(new panzer::IntegrationRule(desc.integrationOrder,m_cell_data));
894 
895  // this function always creates a residual and scatter
896  desc.residualName.first = true;
897 
898  if (residualName == "")
899  desc.residualName.second = "RESIDUAL_" + dofName;
900  else
901  desc.residualName.second = residualName;
902 
903  if (scatterName == "")
904  desc.scatterName = "SCATTER_" + dofName;
905  else
906  desc.scatterName = scatterName;
907 
908 }
909 
910 // ***********************************************************************
911 template <typename EvalT>
913 addDOFGrad(const std::string & dofName,
914  const std::string & gradName)
915 {
916  typename std::map<std::string,DOFDescriptor>::iterator itr = m_provided_dofs_desc.find(dofName);
917 
918  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
919  "EquationSet_DefaultImpl::addDOFGrad: DOF \"" << dofName << "\" has not been specified as a DOF "
920  "by derived equation set \"" << this->getType() << "\".");
921 
922  // allocate and populate a dof descriptor associated with the field "dofName"
923  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
924  TEUCHOS_ASSERT(desc.dofName==dofName); // safety check
925 
926  if (gradName == "")
927  desc.grad = std::make_pair(true,std::string("GRAD_")+dofName);
928  else
929  desc.grad = std::make_pair(true,gradName);
930 }
931 
932 // ***********************************************************************
933 template <typename EvalT>
935 addDOFCurl(const std::string & dofName,
936  const std::string & curlName)
937 {
938  typename std::map<std::string,DOFDescriptor>::iterator itr = m_provided_dofs_desc.find(dofName);
939 
940  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
941  "EquationSet_DefaultImpl::addDOFCurl: DOF \"" << dofName << "\" has not been specified as a DOF "
942  "by derived equation set \"" << this->getType() << "\".");
943 
944  // allocate and populate a dof descriptor associated with the field "dofName"
945  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
946  TEUCHOS_ASSERT(desc.dofName==dofName); // safety check
947 
948  if (curlName == "")
949  desc.curl = std::make_pair(true,std::string("CURL_")+dofName);
950  else
951  desc.curl = std::make_pair(true,curlName);
952 }
953 
954 // ***********************************************************************
955 template <typename EvalT>
957 addDOFDiv(const std::string & dofName,
958  const std::string & divName)
959 {
960  typename std::map<std::string,DOFDescriptor>::iterator itr = m_provided_dofs_desc.find(dofName);
961 
962  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
963  "EquationSet_DefaultImpl::addDOFDiv: DOF \"" << dofName << "\" has not been specified as a DOF "
964  "by derived equation set \"" << this->getType() << "\".");
965 
966  // allocate and populate a dof descriptor associated with the field "dofName"
967  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
968  TEUCHOS_ASSERT(desc.dofName==dofName); // safety check
969 
970  if (divName == "")
971  desc.div = std::make_pair(true,std::string("DIV_")+dofName);
972  else
973  desc.div = std::make_pair(true,divName);
974 }
975 
976 // ***********************************************************************
977 template <typename EvalT>
979 addDOFTimeDerivative(const std::string & dofName,
980  const std::string & dotName)
981 {
982  typename std::map<std::string,DOFDescriptor>::iterator itr = m_provided_dofs_desc.find(dofName);
983 
984  TEUCHOS_TEST_FOR_EXCEPTION(itr==m_provided_dofs_desc.end(),std::runtime_error,
985  "EquationSet_DefaultImpl::addDOFTimeDerivative: DOF \"" << dofName << "\" has not been specified as a DOF "
986  "by derived equation set \"" << this->getType() << "\".");
987 
988  // allocate and populate a dof descriptor associated with the field "dofName"
989  DOFDescriptor & desc = m_provided_dofs_desc[dofName];
990  TEUCHOS_ASSERT(desc.dofName==dofName); // safety check
991 
992  if (dotName == "")
993  desc.timeDerivative = std::make_pair(true,std::string("DXDT_")+dofName);
994  else
995  desc.timeDerivative = std::make_pair(true,dotName);
996 }
997 
998 // ***********************************************************************
999 template <typename EvalT>
1001 setCoordinateDOFs(const std::vector<std::string> & dofNames)
1002 {
1003  TEUCHOS_TEST_FOR_EXCEPTION(m_cell_data.baseCellDimension()!=Teuchos::as<int>(dofNames.size()),std::invalid_argument,
1004  "EquationSet_DefaultImpl::setCoordinateDOFs: Size of vector is not equal to the "
1005  "spatial dimension.");
1006 
1007  for(std::size_t d=0;d<dofNames.size();d++) {
1008  typename std::map<std::string,DOFDescriptor>::const_iterator desc_it = m_provided_dofs_desc.find(dofNames[d]);
1009  TEUCHOS_TEST_FOR_EXCEPTION(desc_it == m_provided_dofs_desc.end(),std::invalid_argument,
1010  "EquationSet_DefaultImpl::setCoordinateDOFs: DOF of name \"" + dofNames[d] + "\" "
1011  "has not been added, thus cannot be set as a coordinate DOF.");
1012  }
1013 
1014  m_coordinate_dofs.push_back(dofNames);
1015 }
1016 
1017 // ***********************************************************************
1018 template <typename EvalT>
1021 {
1022  std::string default_type = "";
1023  valid_parameters.set("Type",default_type,"The equation set type. This must corespond to the type keyword used to build the equation set in the equation set factory.");
1024 }
1025 
1026 // ***********************************************************************
1027 template <typename EvalT>
1030 {
1031  typename std::map<std::string,DOFDescriptor>::const_iterator desc_it = m_provided_dofs_desc.find(dof_name);
1032  TEUCHOS_ASSERT(desc_it != m_provided_dofs_desc.end());
1033  return desc_it->second.basis;
1034 }
1035 
1036 // ***********************************************************************
1037 template <typename EvalT>
1040 {
1041  typename std::map<std::string,DOFDescriptor>::const_iterator desc_it = m_provided_dofs_desc.find(dof_name);
1042  TEUCHOS_TEST_FOR_EXCEPTION(desc_it == m_provided_dofs_desc.end(),std::logic_error,
1043  "EquationSet_DefaultImpl::getIntRuleForDOF: Failed to find degree of freedom "
1044  "with name \"" << dof_name << "\".");
1045  return desc_it->second.intRule;
1046 }
1047 
1048 // ***********************************************************************
1049 template <typename EvalT>
1052 {
1053  typename std::map<std::string,DOFDescriptor>::const_iterator desc_it = m_provided_dofs_desc.find(dof_name);
1054  TEUCHOS_TEST_FOR_EXCEPTION(desc_it == m_provided_dofs_desc.end(),std::logic_error,
1055  "EquationSet_DefaultImpl::getBasisIRLayoutForDOF: Failed to find degree of freedom "
1056  "with name \"" << dof_name << "\".");
1057 
1058  return panzer::basisIRLayout(desc_it->second.basis,*(desc_it->second.intRule));
1059 }
1060 
1061 // ***********************************************************************
1062 template <typename EvalT>
1065  const std::string dof_name,
1066  const std::vector<std::string>& residual_contributions,
1067  const std::string residual_field_name) const
1068 {
1069  using Teuchos::rcp;
1070  using Teuchos::RCP;
1071 
1073 
1074  if (residual_field_name != "")
1075  p.set("Sum Name", residual_field_name);
1076  else
1077  p.set("Sum Name", "RESIDUAL_" + dof_name);
1078 
1079  RCP<std::vector<std::string> > rcp_residual_contributions = rcp(new std::vector<std::string>);
1080  *rcp_residual_contributions = residual_contributions;
1081 
1082  p.set("Values Names", rcp_residual_contributions);
1083 
1084  DescriptorIterator desc_it = m_provided_dofs_desc.find(dof_name);
1085  TEUCHOS_ASSERT(desc_it != m_provided_dofs_desc.end());
1086 
1087  p.set("Data Layout", desc_it->second.basis->functional);
1088 
1090 
1091  this->template registerEvaluator<EvalT>(fm, op);
1092 }
1093 
1094 // ***********************************************************************
1095 template <typename EvalT>
1098  const std::string dof_name,
1099  const std::vector<std::string>& residual_contributions,
1100  const std::vector<double>& scale_contributions,
1101  const std::string residual_field_name) const
1102 {
1103  using Teuchos::rcp;
1104  using Teuchos::RCP;
1105 
1107 
1108  if (residual_field_name != "")
1109  p.set("Sum Name", residual_field_name);
1110  else
1111  p.set("Sum Name", "RESIDUAL_" + dof_name);
1112 
1113  RCP<std::vector<std::string> > rcp_residual_contributions = rcp(new std::vector<std::string>);
1114  *rcp_residual_contributions = residual_contributions;
1115  p.set("Values Names", rcp_residual_contributions);
1116 
1117  RCP<const std::vector<double> > rcp_scale_contributions = rcp(new std::vector<double>(scale_contributions));
1118  p.set("Scalars", rcp_scale_contributions);
1119 
1120  DescriptorIterator desc_it = m_provided_dofs_desc.find(dof_name);
1121  TEUCHOS_ASSERT(desc_it != m_provided_dofs_desc.end());
1122 
1123  p.set("Data Layout", desc_it->second.basis->functional);
1124 
1126 
1127  this->template registerEvaluator<EvalT>(fm, op);
1128 }
1129 
1130 // ***********************************************************************
1131 template <typename EvalT>
1132 void panzer::EquationSet_DefaultImpl<EvalT>::addClosureModel(const std::string& closure_model)
1133 {
1134  m_closure_model_ids.push_back(closure_model);
1135 }
1136 
1137 // ***********************************************************************
1138 template <typename EvalT>
1141 {
1142  return m_input_params;
1143 }
1144 
1145 // ***********************************************************************
1146 #endif
void buildAndRegisterResidualSummationEvaluator(PHX::FieldManager< panzer::Traits > &fm, const std::string dof_name, const std::vector< std::string > &residual_contributions, const std::string residual_field_name="") const
Teuchos::RCP< panzer::BasisIRLayout > lookupLayout(const std::string &fieldName) const
Get the basis associated with a particular field.
Teuchos::RCP< panzer::IntegrationRule > getIntRuleForDOF(const std::string &dof_name) const
Returns the integration rule associated with the residual contributions for the dof_name.
std::map< std::string, DOFDescriptor >::const_iterator DescriptorIterator
For convenience, declare the DOFDescriptor iterator.
Interpolates basis DOF values to IP DOF Gradient values.
void addDOFGrad(const std::string &dofName, const std::string &gradName="")
bool is_null(const boost::shared_ptr< T > &p)
virtual void setTangentParamNames(const std::vector< std::string > &tangent_param_names)
Set the list of tangent parameter names.
virtual std::string getType() const
Returns the type of the equation set object. Corresponds to the keyword used by the equation set fact...
Teuchos::RCP< typename Sacado::mpl::apply< panzer::ClosureModelFactory< _ >, ScalarT >::type > getAsObject()
int getIntegrationOrder(const std::string &dofName) const
Get the integration order for an existing degree of freedom.
Default implementation for accessing the GlobalData object.
virtual const std::map< int, Teuchos::RCP< panzer::IntegrationRule > > & getIntegrationRules() const
Return a map of unique integration rules for the equation set, key is the integration order...
T & get(const std::string &name, T def_value)
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 void buildAndRegisterInitialConditionEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const std::string &model_name, const Teuchos::ParameterList &models, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
Teuchos::RCP< Teuchos::ParameterList > m_eval_plist
virtual void buildAndRegisterDOFProjectionsToIPEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLayoutLibrary &fl, const Teuchos::RCP< panzer::IntegrationRule > &ir, const Teuchos::Ptr< const panzer::LinearObjFactory< panzer::Traits > > &lof, const Teuchos::ParameterList &user_data) const
int getBasisOrder(const std::string &dofName) const
Get the basis order for an existing degree of freedom.
Teuchos::RCP< PHX::Evaluator< Traits > > buildGatherTangent(const Teuchos::ParameterList &pl) const
Use preconstructed gather evaluators.
Teuchos::RCP< panzer::PureBasis > getBasisForDOF(const std::string &dof_name) const
Returns the PureBasis associated with the residual contributions for the dof_name.
Interpolates basis DOF values to IP DOF Div values.
void setDefaultValidParameters(Teuchos::ParameterList &valid_parameters)
PHX::View< const int * > offsets
Teuchos::RCP< Teuchos::ParameterList > getEquationSetParameterList() const
Returns the parameter list used to build this equation set.
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
bool isParameter(const std::string &name) const
void setCoordinateDOFs(const std::vector< std::string > &dofNames)
Teuchos::RCP< PHX::Evaluator< Traits > > buildScatterDirichlet(const Teuchos::ParameterList &pl) const
Use preconstructed dirichlet scatter evaluators.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Data for determining cell topology and dimensionality.
void addDOFDiv(const std::string &dofName, const std::string &divName="")
virtual const Teuchos::RCP< Teuchos::ParameterList > getEvaluatorParameterList() const
Returns the parameter list that will be passed off from the equaiton set to the closure model evaluat...
Gathers coordinates for the quadrature from the workset and stores them in the field manager...
virtual void setupDOFs()
Builds the integration rule, basis, DOFs, and default parameter list. This MUST be called in the cons...
std::map< std::string, std::pair< Teuchos::RCP< panzer::PureBasis >, Teuchos::RCP< std::vector< std::string > > > >::const_iterator BasisIterator
For convenience, declare a basis iterator.
Interpolates basis DOF values to IP DOF Curl values.
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
Teuchos::RCP< PHX::Evaluator< Traits > > buildGather(const Teuchos::ParameterList &pl) const
Use preconstructed gather evaluators.
void addDOF(const std::string &dofName, const std::string &basisType, const int &basisOrder, const int integrationOrder=-1, const std::string residualName="", const std::string scatterName="")
void updateDOF(const std::string &dofName, int basisOrder, int integrationOrder=-1)
Modifying an existing DOF&#39;s basis function and integration rule.
virtual void buildAndRegisterClosureModelEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLayoutLibrary &fl, const Teuchos::RCP< panzer::IntegrationRule > &ir, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const Teuchos::ParameterList &models, const Teuchos::ParameterList &user_data) const
Register closure model evaluators with the model name internally specified by the equation set...
Teuchos::RCP< panzer::BasisIRLayout > getBasisIRLayoutForDOF(const std::string &dof_name) const
Returns the BasisIRLayout for the dof_name.
virtual void buildAndRegisterGatherAndOrientationEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
bool nonnull(const boost::shared_ptr< T > &p)
Teuchos::RCP< const FieldLayoutLibrary > buildFieldLayoutLibrary(panzer::PointRule &ir) const
Teuchos::RCP< PHX::Evaluator< Traits > > buildGatherOrientation(const Teuchos::ParameterList &pl) const
Use preconstructed gather evaluators.
const Teuchos::RCP< Teuchos::ParameterList > m_input_params
virtual const std::vector< int > & getGIDFieldOffsets(const std::string &blockId, int fieldNum) const =0
Use the field pattern so that you can find a particular field in the GIDs array.
void getAddedDOFs(std::vector< std::string > &dofNames) const
virtual const std::vector< std::pair< std::string, Teuchos::RCP< panzer::PureBasis > > > & getProvidedDOFs() const
Return the Basis for the equation set, key is the DOF name (note coordinate DOFs are NOT included) ...
Interpolates basis DOF values to IP DOF values.
Definition: Panzer_DOF.hpp:56
virtual const std::vector< std::vector< std::string > > & getCoordinateDOFs() const
Return a vector of vectors that correspond to DOFs set as coordinate fields.
void addDOFTimeDerivative(const std::string &dofName, const std::string &dotName="")
Teuchos::RCP< PHX::Evaluator< Traits > > buildScatter(const Teuchos::ParameterList &pl) const
Use preconstructed scatter evaluators.
bool buildTransientSupport() const
Returns true if transient support should be enabled in the equation set.
virtual void buildAndRegisterScatterEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
void addDOFCurl(const std::string &dofName, const std::string &curlName="")
EquationSet_DefaultImpl(const Teuchos::RCP< Teuchos::ParameterList > &params, const int &default_integration_order, const panzer::CellData &cell_data, const Teuchos::RCP< panzer::GlobalData > &global_data, const bool build_transient_support)
virtual void setElementBlockId(const std::string &blockId)
Gathers coordinates for the basis function from the workset and stores them in the field manager...
void addClosureModel(const std::string &closure_model)