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 //
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_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
45 
46 #include "Intrepid2_FunctionSpaceTools.hpp"
48 #include "Panzer_BasisIRLayout.hpp"
50 #include "Kokkos_ViewFactory.hpp"
51 
52 namespace panzer {
53 
54 //**********************************************************************
55 template<typename EvalT, typename Traits>
58  const Teuchos::ParameterList& p) :
59  residual( p.get<std::string>("Residual Name"),
60  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
61  scalar( p.get<std::string>("Value Name"),
62  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar),
63  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
64 {
66  p.validateParameters(*valid_params);
67 
69  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
70 
71  // Verify that this basis supports the gradient operation
72  TEUCHOS_TEST_FOR_EXCEPTION(!basis->isScalarBasis(),std::logic_error,
73  "Integrator_TransientBasisTimesScalar: Basis of type \"" << basis->name() << "\" is not a "
74  "scalar basis");
75 
76  this->addEvaluatedField(residual);
77  this->addDependentField(scalar);
78 
79  multiplier = p.get<double>("Multiplier");
80 
81 
82  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
83  const std::vector<std::string>& field_multiplier_names =
84  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
85 
86  for (std::vector<std::string>::const_iterator name =
87  field_multiplier_names.begin();
88  name != field_multiplier_names.end(); ++name) {
89  PHX::MDField<const ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
90  field_multipliers.push_back(tmp_field);
91  }
92  }
93 
94  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
95  field != field_multipliers.end(); ++field)
96  this->addDependentField(*field);
97 
98  std::string n = "Integrator_TransientBasisTimesScalar: " + residual.fieldTag().name();
99  this->setName(n);
100 }
101 
102 //**********************************************************************
103 template<typename EvalT, typename Traits>
104 void
107  typename Traits::SetupData sd,
109 {
110  this->utils.setFieldData(residual,fm);
111  this->utils.setFieldData(scalar,fm);
112 
113  num_nodes = residual.extent(1);
114  num_qp = scalar.extent(1);
115 
116  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
117 
118  tmp = Kokkos::createDynRankView(residual.get_static_view(),"tmp",scalar.extent(0), num_qp);
119 }
120 
121 //**********************************************************************
122 template<typename EvalT, typename Traits>
123 void
126  typename Traits::EvalData workset)
127 {
128  if (workset.evaluate_transient_terms) {
129 
130  // for (int i=0; i < residual.size(); ++i)
131  // residual[i] = 0.0;
132 
133  Kokkos::deep_copy (residual.get_static_view(), ScalarT(0.0));
134 
135  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
136  for (std::size_t qp = 0; qp < num_qp; ++qp) {
137  tmp(cell,qp) = multiplier * scalar(cell,qp);
138  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
139  field != field_multipliers.end(); ++field)
140  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
141  }
142  }
143 
144  if(workset.num_cells>0)
145  Intrepid2::FunctionSpaceTools<PHX::exec_space>::
146  integrate<ScalarT>(residual.get_view(),
147  tmp,
148  (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
149  }
150 }
151 
152 //**********************************************************************
153 template<typename EvalT, typename TRAITS>
156 {
158  p->set<std::string>("Residual Name", "?");
159  p->set<std::string>("Value Name", "?");
161  p->set("Basis", basis);
163  p->set("IR", ir);
164  p->set<double>("Multiplier", 1.0);
166  p->set("Field Multipliers", fms);
167  return p;
168 }
169 
170 //**********************************************************************
171 
172 }
173 
174 #endif
175 
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.
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)
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)
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_