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<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_.dimension(1);
67  int numPoints = grad_basis_.dimension(2);
68  int spaceDim = grad_basis_.dimension(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<ScalarT,Cell,Point> dof_value_;
81  const ArrayT grad_basis_;
82 
83 };
84 
85 }
86 
87 //**********************************************************************
88 PHX_EVALUATOR_CTOR(DOFGradient,p) :
89  dof_value( p.get<std::string>("Name"),
90  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
91  dof_gradient( p.get<std::string>("Gradient Name"),
92  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_vector ),
93  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
94 {
96  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
97 
98  // Verify that this basis supports the gradient operation
99  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
100  "DOFGradient: Basis of type \"" << basis->name() << "\" does not support GRAD");
101 
102  this->addEvaluatedField(dof_gradient);
103  this->addDependentField(dof_value);
104 
105  std::string n = "DOFGradient: " + dof_gradient.fieldTag().name() + " ("+PHX::typeAsString<EvalT>()+")";
106  this->setName(n);
107 }
108 
109 //**********************************************************************
110 PHX_POST_REGISTRATION_SETUP(DOFGradient,sd,fm)
111 {
112  this->utils.setFieldData(dof_value,fm);
113  this->utils.setFieldData(dof_gradient,fm);
114 
115  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
116 }
117 
118 //**********************************************************************
119 PHX_EVALUATE_FIELDS(DOFGradient,workset)
120 {
121 /*
122  // Zero out arrays (probably don't need this anymore)
123  dof_gradient.deep_copy(ScalarT(0.0));
124 
125  if(workset.num_cells>0)
126  Intrepid2::FunctionSpaceTools::evaluate<ScalarT>(dof_gradient,dof_value,(this->wda(workset).bases[basis_index])->grad_basis);
127 */
128  if (workset.num_cells == 0 )
129  return;
130  typedef decltype(this->wda(workset).bases[basis_index]->grad_basis) ArrayT;
131  evaluateGrad_withSens<ScalarT, ArrayT> eval(dof_gradient,dof_value,this->wda(workset).bases[basis_index]->grad_basis);
132  Kokkos::parallel_for(workset.num_cells, eval);
133 }
134 
135 //**********************************************************************
136 
137 }
138 
139 #endif
std::size_t basis_index
const ArrayT grad_basis_
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, Point > dof_value_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void operator()(const FieldMultTag< NUM_FIELD_MULT > &, const std::size_t &cell) const
PHX::MDField< ScalarT > dof_grad_
bool supportsGrad() const
int numPoints
T * get() const
int numFields
PHX::MDField< ScalarT > dof_gradient
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
PHX::MDField< ScalarT, Cell, Point > dof_value
Interpolates basis DOF values to IP DOF Gradient values.
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)