Panzer  Version of the Day
Panzer_Integrator_DivBasisTimesScalar_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_Integrator_DivBasisTimesScalar_impl_hpp__
44 #define __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__
45 
46 #include "Intrepid2_FunctionSpaceTools.hpp"
47 
49 #include "Panzer_BasisIRLayout.hpp"
51 #include "Kokkos_ViewFactory.hpp"
52 
53 namespace panzer {
54 
55 //**********************************************************************
56 template<typename EvalT, typename Traits>
59  const Teuchos::ParameterList& p) :
60  residual( p.get<std::string>("Residual Name"),
61  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
62  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
63 {
64  Teuchos::RCP<Teuchos::ParameterList> valid_params = this->getValidParameters();
65  p.validateParameters(*valid_params);
66 
67  Teuchos::RCP<const PureBasis> basis
68  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
69 
70  // Verify that this basis supports the curl operation
71  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsDiv(),std::logic_error,
72  "Integrator_DivBasisTimesScalar: Basis of type \"" << basis->name() << "\" does not support DIV.");
73  TEUCHOS_TEST_FOR_EXCEPTION(!basis->requiresOrientations(),std::logic_error,
74  "Integration_DivBasisTimesScalar: Basis of type \"" << basis->name() << "\" should require orientations. So we are throwing.");
75 
76  scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Value Name"),
77  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar );
78 
79  this->addEvaluatedField(residual);
80  this->addDependentField(scalar);
81 
82  multiplier = p.get<double>("Multiplier");
83  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
84 
85  const std::vector<std::string>& field_multiplier_names =
86  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
87 
88  for (const std::string & name : field_multiplier_names) {
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 (const auto & field : field_multipliers)
95  this->addDependentField(field);
96 
97  std::string n =
98  "Integrator_DivBasisTimesScalar: " + residual.fieldTag().name();
99 
100  this->setName(n);
101 }
102 
103 //**********************************************************************
104 template<typename EvalT, typename Traits>
105 void
108  typename Traits::SetupData sd,
110 {
111  this->utils.setFieldData(residual,fm);
112  this->utils.setFieldData(scalar,fm);
113 
114  for (auto & field : field_multipliers)
115  this->utils.setFieldData(field,fm);
116 
117  num_nodes = residual.extent(1);
118  num_qp = scalar.extent(1);
119 
120  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
121 
122  tmp = Kokkos::createDynRankView(residual.get_static_view(),"tmp",scalar.extent(0), num_qp);
123 }
124 
125 //**********************************************************************
126 template<typename EvalT, typename Traits>
127 void
130  typename Traits::EvalData workset)
131 {
132  // zero the reisdual
133  residual.deep_copy(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  ScalarT tmpVar = 1.0;
138  for (const auto & field : field_multipliers)
139  tmpVar = tmpVar * field(cell,qp);
140 
141  // no dimension to loop over for scalar fields
142  tmp(cell,qp) = multiplier * tmpVar * scalar(cell,qp);
143  }
144  }
145 
146  {
147  // const Kokkos::DynRankView<double,PHX::Device> & weighted_div_basis = (this->wda(workset).bases[basis_index])->weighted_div_basis;
148  const BasisValues2<double> & bv = *this->wda(workset).bases[basis_index];
149 
150  for (index_t cell = 0; cell < workset.num_cells; ++cell)
151  for (std::size_t basis = 0; basis < num_nodes; ++basis) {
152  for (std::size_t qp = 0; qp < num_qp; ++qp)
153  residual(cell,basis) += tmp(cell,qp)*bv.weighted_div_basis(cell,basis,qp);
154  }
155  }
156 /*
157  if(workset.num_cells>0) {
158  Intrepid2::FunctionSpaceTools::
159  integrate<ScalarT>(residual, tmp,
160  this->wda(workset).bases[basis_index]->weighted_div_basis,
161  Intrepid2::COMP_CPP);
162  }
163 */
164 }
165 
166 //**********************************************************************
167 
168 template<typename EvalT, typename TRAITS>
169 Teuchos::RCP<Teuchos::ParameterList>
171 {
172  Teuchos::RCP<Teuchos::ParameterList> p = Teuchos::rcp(new Teuchos::ParameterList);
173  p->set<std::string>("Residual Name", "?");
174  p->set<std::string>("Value Name", "?");
175  p->set<std::string>("Test Field Name", "?");
176  Teuchos::RCP<panzer::BasisIRLayout> basis;
177  p->set("Basis", basis);
178  Teuchos::RCP<panzer::IntegrationRule> ir;
179  p->set("IR", ir);
180  p->set<double>("Multiplier", 1.0);
181  Teuchos::RCP<const std::vector<std::string> > fms;
182  p->set("Field Multipliers", fms);
183  return p;
184 }
185 
186 //**********************************************************************
187 
188 }
189 
190 #endif
double multiplier
The scalar multiplier out in front of the integral ( ).
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.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
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.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
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...
Array_CellBasisIP weighted_div_basis
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_