Panzer  Version of the Day
Panzer_DOFGradient_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_DOF_GRADIENT_IMPL_HPP
44 #define PANZER_DOF_GRADIENT_IMPL_HPP
45 
47 #include "Panzer_BasisIRLayout.hpp"
49 #include "Intrepid2_FunctionSpaceTools.hpp"
50 
51 namespace panzer {
52 
53 namespace {
54 
55 template <typename ScalarT,typename ArrayT>
56 struct evaluateGrad_withSens {
57  evaluateGrad_withSens (PHX::MDField<ScalarT> dof_grad,
58  PHX::MDField<const ScalarT,Cell,Point> dof_value,
59  const ArrayT grad_basis) :
60  dof_grad_(dof_grad),dof_value_(dof_value),grad_basis_(grad_basis)
61  {}
62  KOKKOS_INLINE_FUNCTION
63  void operator() (const size_t &cell) const
64  {
65  // evaluate at quadrature points
66  int numFields = grad_basis_.extent(1);
67  int numPoints = grad_basis_.extent(2);
68  int spaceDim = grad_basis_.extent(3);
69  for (int pt=0; pt<numPoints; pt++) {
70  for (int d=0; d<spaceDim; d++) {
71  // first initialize to the right thing (prevents over writing with 0)
72  // then loop over one less basis function
73  dof_grad_(cell,pt,d) = dof_value_(cell, 0) * grad_basis_(cell, 0, pt, d);
74  for (int bf=1; bf<numFields; bf++)
75  dof_grad_(cell,pt,d) += dof_value_(cell, bf) * grad_basis_(cell, bf, pt, d);
76  }
77  }
78  }
79  PHX::MDField<ScalarT> dof_grad_;
80  PHX::MDField<const ScalarT,Cell,Point> dof_value_;
81  const ArrayT grad_basis_;
82 
83 };
84 
85 }
86 
87 //**********************************************************************
88 template<typename EvalT, typename Traits>
91  const Teuchos::ParameterList& p) :
92  use_descriptors_(false),
93  dof_value( p.get<std::string>("Name"),
94  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
95  dof_gradient( p.get<std::string>("Gradient Name"),
96  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_vector ),
97  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
98 {
100  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
101 
102  // Verify that this basis supports the gradient operation
103  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
104  "DOFGradient: Basis of type \"" << basis->name() << "\" does not support GRAD");
105 
106  this->addEvaluatedField(dof_gradient);
107  this->addDependentField(dof_value);
108 
109  std::string n = "DOFGradient: " + dof_gradient.fieldTag().name() + " ("+PHX::typeAsString<EvalT>()+")";
110  this->setName(n);
111 }
112 
113 //**********************************************************************
114 template<typename EvalT, typename TRAITS>
116 DOFGradient(const PHX::FieldTag & input,
117  const PHX::FieldTag & output,
118  const panzer::BasisDescriptor & bd,
120  : use_descriptors_(true)
121  , bd_(bd)
122  , id_(id)
123  , dof_value(input)
124  , dof_gradient(output)
125 {
126  // Verify that this basis supports the gradient operation
127  TEUCHOS_TEST_FOR_EXCEPTION(bd_.getType()=="HGrad",std::logic_error,
128  "DOFGradient: Basis of type \"" << bd_.getType() << "\" does not support GRAD");
129 
130  this->addEvaluatedField(dof_gradient);
131  this->addDependentField(dof_value);
132 
133  std::string n = "DOFGradient: " + dof_gradient.fieldTag().name() + " ("+PHX::typeAsString<EvalT>()+")";
134  this->setName(n);
135 }
136 
137 //**********************************************************************
138 template<typename EvalT, typename Traits>
139 void
142  typename Traits::SetupData sd,
144 {
145  this->utils.setFieldData(dof_value,fm);
146  this->utils.setFieldData(dof_gradient,fm);
147 
148  if(not use_descriptors_)
149  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
150 }
151 
152 //**********************************************************************
153 template<typename EvalT, typename Traits>
154 void
157  typename Traits::EvalData workset)
158 {
159  if (workset.num_cells == 0 )
160  return;
161 
162  const panzer::BasisValues2<double> & basisValues = use_descriptors_ ? this->wda(workset).getBasisValues(bd_,id_)
163  : *this->wda(workset).bases[basis_index];
164 
165  typedef decltype(basisValues.grad_basis) ArrayT;
166 
167  evaluateGrad_withSens<ScalarT, ArrayT> eval(dof_gradient,dof_value,basisValues.grad_basis);
168  Kokkos::parallel_for(workset.num_cells, eval);
169 }
170 
171 //**********************************************************************
172 
173 }
174 
175 #endif
const ArrayT grad_basis_
T & get(const std::string &name, T def_value)
DOFGradient(const Teuchos::ParameterList &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
PHX::MDField< ScalarT > dof_gradient
PHX::MDField< ScalarT > dof_grad_
int numPoints
int numFields
panzer::BasisDescriptor bd_
PHX::MDField< const ScalarT, Cell, Point > dof_value
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 evaluateFields(typename TRAITS::EvalData d)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
PHX::MDField< const ScalarT, Cell, Point > dof_value_
Array_CellBasisIPDim grad_basis
PHX::MDField< const ScalarT, Cell, Point > dof_value
const std::string & getType() const
Get type of basis.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_