Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_Scalar_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_SCALAR_IMPL_HPP
44#define PANZER_EVALUATOR_SCALAR_IMPL_HPP
45
46#include "Intrepid2_FunctionSpaceTools.hpp"
49#include "Kokkos_ViewFactory.hpp"
50#include "Phalanx_DataLayout_MDALayout.hpp"
51
52namespace panzer {
53
54//**********************************************************************
55template<typename EvalT, typename Traits>
58 const Teuchos::ParameterList& p) : quad_index(0)
59{
60 Teuchos::RCP<Teuchos::ParameterList> valid_params = this->getValidParameters();
61 p.validateParameters(*valid_params);
62
63 Teuchos::RCP<panzer::IntegrationRule> ir = p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
64 quad_order = ir->cubature_degree;
65
66 Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(ir->dl_scalar->extent(0)));
67 integral = PHX::MDField<ScalarT>( p.get<std::string>("Integral Name"), dl_cell);
68 scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Integrand Name"), ir->dl_scalar);
69
70 this->addEvaluatedField(integral);
71 this->addDependentField(scalar);
72
73 multiplier = 1.0;
74 if(p.isType<double>("Multiplier"))
75 multiplier = p.get<double>("Multiplier");
76
77 if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
78 const std::vector<std::string>& field_multiplier_names =
79 *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
80
81 field_multipliers_h = typename PHX::View<PHX::UnmanagedView<const ScalarT**>* >::HostMirror("FieldMultipliersHost", field_multiplier_names.size());
82 field_multipliers = PHX::View<PHX::UnmanagedView<const ScalarT**>* >("FieldMultipliersDevice", field_multiplier_names.size());
83
84 int cnt=0;
85 for (std::vector<std::string>::const_iterator name =
86 field_multiplier_names.begin();
87 name != field_multiplier_names.end(); ++name, ++cnt) {
88 PHX::MDField<const ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
89 this->addDependentField(tmp_field.fieldTag(),field_multipliers_h(cnt));
90 }
91 }
92
93 std::string n = "Integrator_Scalar: " + integral.fieldTag().name();
94 this->setName(n);
95}
96
97//**********************************************************************
98template<typename EvalT, typename Traits>
99void
102 typename Traits::SetupData sd,
104{
105 num_qp = scalar.extent(1);
106 tmp = Kokkos::createDynRankView(scalar.get_static_view(),"tmp", scalar.extent(0), num_qp);
107 quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
108
109 Kokkos::deep_copy(field_multipliers, field_multipliers_h);
110}
111
112
113//**********************************************************************
114template<typename EvalT, typename Traits>
115void
118 typename Traits::EvalData workset)
119{
120
121 const auto wm = this->wda(workset).int_rules[quad_index]->weighted_measure;
122
123 auto l_tmp = tmp;
124 auto l_multiplier = multiplier;
125 auto l_scalar = scalar;
126 auto l_field_multipliers = field_multipliers;
127 auto l_num_qp = num_qp;
128 auto l_integral = integral.get_static_view();
129
130 Kokkos::parallel_for("IntegratorScalar", workset.num_cells, KOKKOS_LAMBDA (int cell) {
131 for (std::size_t qp = 0; qp < l_num_qp; ++qp) {
132 l_tmp(cell,qp) = l_multiplier * l_scalar(cell,qp);
133 for (size_t f=0;f<l_field_multipliers.extent(0); ++f)
134 l_tmp(cell,qp) *= l_field_multipliers(f)(cell,qp);
135 }
136 l_integral(cell) = 0.0;
137 for (std::size_t qp = 0; qp < l_num_qp; ++qp) {
138 l_integral(cell) += l_tmp(cell, qp)*wm(cell, qp);
139 }
140 } );
141 Kokkos::fence();
142}
143
144//**********************************************************************
145template<typename EvalT, typename TRAITS>
146Teuchos::RCP<Teuchos::ParameterList>
148{
149 Teuchos::RCP<Teuchos::ParameterList> p = Teuchos::rcp(new Teuchos::ParameterList);
150 p->set<std::string>("Integral Name", "?");
151 p->set<std::string>("Integrand Name", "?");
152
153 Teuchos::RCP<panzer::IntegrationRule> ir;
154 p->set("IR", ir);
155 p->set<double>("Multiplier", 1.0);
156
157 Teuchos::RCP<const std::vector<std::string> > fms;
158 p->set("Field Multipliers", fms);
159 return p;
160}
161
162//**********************************************************************
163
164}
165
166#endif
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< const ScalarT, Cell, IP > scalar
Integrator_Scalar(const Teuchos::ParameterList &p)
PHX::View< PHX::UnmanagedView< constScalarT ** > * >::HostMirror field_multipliers_h
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > field_multipliers
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
int num_cells
DEPRECATED - use: numCells()
std::vector< int >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_