Panzer  Version of the Day
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 PHX_EVALUATOR_CTOR(Integrator_TransientBasisTimesScalar,p) :
56  residual( p.get<std::string>("Residual Name"),
57  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
58  scalar( p.get<std::string>("Value Name"),
59  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar),
60  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
61 {
63  p.validateParameters(*valid_params);
64 
66  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
67 
68  // Verify that this basis supports the gradient operation
69  TEUCHOS_TEST_FOR_EXCEPTION(!basis->isScalarBasis(),std::logic_error,
70  "Integrator_TransientBasisTimesScalar: Basis of type \"" << basis->name() << "\" is not a "
71  "scalar basis");
72 
73  this->addEvaluatedField(residual);
74  this->addDependentField(scalar);
75 
76  multiplier = p.get<double>("Multiplier");
77 
78 
79  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
80  const std::vector<std::string>& field_multiplier_names =
81  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
82 
83  for (std::vector<std::string>::const_iterator name =
84  field_multiplier_names.begin();
85  name != field_multiplier_names.end(); ++name) {
86  PHX::MDField<ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
87  field_multipliers.push_back(tmp_field);
88  }
89  }
90 
91  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
92  field != field_multipliers.end(); ++field)
93  this->addDependentField(*field);
94 
95  std::string n = "Integrator_TransientBasisTimesScalar: " + residual.fieldTag().name();
96  this->setName(n);
97 }
98 
99 //**********************************************************************
100 PHX_POST_REGISTRATION_SETUP(Integrator_TransientBasisTimesScalar,sd,fm)
101 {
102  this->utils.setFieldData(residual,fm);
103  this->utils.setFieldData(scalar,fm);
104 
105  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
106  field != field_multipliers.end(); ++field)
107  this->utils.setFieldData(*field,fm);
108 
109  num_nodes = residual.dimension(1);
110  num_qp = scalar.dimension(1);
111 
112  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
113 
114  tmp = Kokkos::createDynRankView(residual.get_static_view(),"tmp",scalar.dimension(0), num_qp);
115 }
116 
117 //**********************************************************************
118 PHX_EVALUATE_FIELDS(Integrator_TransientBasisTimesScalar,workset)
119 {
120  if (workset.evaluate_transient_terms) {
121 
122  // for (int i=0; i < residual.size(); ++i)
123  // residual[i] = 0.0;
124 
125  Kokkos::deep_copy (residual.get_static_view(), ScalarT(0.0));
126 
127  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
128  for (std::size_t qp = 0; qp < num_qp; ++qp) {
129  tmp(cell,qp) = multiplier * scalar(cell,qp);
130  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
131  field != field_multipliers.end(); ++field)
132  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
133  }
134  }
135 
136  if(workset.num_cells>0)
137  Intrepid2::FunctionSpaceTools::
138  integrate<ScalarT>(residual, tmp,
139  (this->wda(workset).bases[basis_index])->weighted_basis_scalar,
140  Intrepid2::COMP_CPP);
141  }
142 }
143 
144 //**********************************************************************
145 template<typename EvalT, typename TRAITS>
148 {
150  p->set<std::string>("Residual Name", "?");
151  p->set<std::string>("Value Name", "?");
153  p->set("Basis", basis);
155  p->set("IR", ir);
156  p->set<double>("Multiplier", 1.0);
158  p->set("Field Multipliers", fms);
159  return p;
160 }
161 
162 //**********************************************************************
163 
164 }
165 
166 #endif
167 
PHX::MDField< ScalarT > residual
Evaluates a Dirichlet BC residual corresponding to a field value.
std::size_t basis_index
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
std::string name() const
A unique key that is the combination of the basis type and basis order.
std::string basis_name
PHX::MDField< ScalarT, Cell > tmp
std::size_t num_qp
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)
T * get() const
PHX::MDField< ScalarT > vector
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.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
bool isScalarBasis() const
double multiplier
PHX::MDField< const ScalarT, Cell, IP > scalar
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)