Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_TransientBasisTimesScalar_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_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
12 #define PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
13 
14 #include "Intrepid2_FunctionSpaceTools.hpp"
16 #include "Panzer_BasisIRLayout.hpp"
18 #include "Kokkos_ViewFactory.hpp"
19 
20 namespace panzer {
21 
22 //**********************************************************************
23 template<typename EvalT, typename Traits>
26  const Teuchos::ParameterList& p) :
27  residual( p.get<std::string>("Residual Name"),
28  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
29  scalar( p.get<std::string>("Value Name"),
30  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar),
31  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
32 {
34  p.validateParameters(*valid_params);
35 
37  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
38 
39  // Verify that this basis supports the gradient operation
40  TEUCHOS_TEST_FOR_EXCEPTION(!basis->isScalarBasis(),std::logic_error,
41  "Integrator_TransientBasisTimesScalar: Basis of type \"" << basis->name() << "\" is not a "
42  "scalar basis");
43 
44  this->addEvaluatedField(residual);
45  this->addDependentField(scalar);
46 
47  multiplier = p.get<double>("Multiplier");
48 
49 
50  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
51  const std::vector<std::string>& field_multiplier_names =
52  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
53 
54  for (std::vector<std::string>::const_iterator name =
55  field_multiplier_names.begin();
56  name != field_multiplier_names.end(); ++name) {
58  field_multipliers.push_back(tmp_field);
59  }
60  }
61 
62  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
63  field != field_multipliers.end(); ++field)
64  this->addDependentField(*field);
65 
66  std::string n = "Integrator_TransientBasisTimesScalar: " + residual.fieldTag().name();
67  this->setName(n);
68 }
69 
70 //**********************************************************************
71 template<typename EvalT, typename Traits>
72 void
75  typename Traits::SetupData sd,
77 {
78  this->utils.setFieldData(residual,fm);
79  this->utils.setFieldData(scalar,fm);
80 
81  num_nodes = residual.extent(1);
82  num_qp = scalar.extent(1);
83 
84  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
85 
86  tmp = Kokkos::createDynRankView(residual.get_static_view(),"tmp",scalar.extent(0), num_qp);
87 }
88 
89 //**********************************************************************
90 template<typename EvalT, typename Traits>
91 void
94  typename Traits::EvalData workset)
95 {
96  if (workset.evaluate_transient_terms) {
97 
98  // for (int i=0; i < residual.size(); ++i)
99  // residual[i] = 0.0;
100 
101  Kokkos::deep_copy (residual.get_static_view(), ScalarT(0.0));
102 
103  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
104  for (std::size_t qp = 0; qp < num_qp; ++qp) {
105  tmp(cell,qp) = multiplier * scalar(cell,qp);
106  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
107  field != field_multipliers.end(); ++field)
108  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
109  }
110  }
111 
112  if(workset.num_cells>0)
113  Intrepid2::FunctionSpaceTools<PHX::exec_space>::
114  integrate<ScalarT>(residual.get_view(),
115  tmp,
116  (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
117  }
118 }
119 
120 //**********************************************************************
121 template<typename EvalT, typename TRAITS>
124 {
126  p->set<std::string>("Residual Name", "?");
127  p->set<std::string>("Value Name", "?");
129  p->set("Basis", basis);
131  p->set("IR", ir);
132  p->set<double>("Multiplier", 1.0);
134  p->set("Field Multipliers", fms);
135  return p;
136 }
137 
138 //**********************************************************************
139 
140 }
141 
142 #endif
143 
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
bool isType(const std::string &name) const
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_