Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_TransientBasisTimesScalar_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
12#define PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
13
14#include "Intrepid2_FunctionSpaceTools.hpp"
18#include "Kokkos_ViewFactory.hpp"
19
20namespace panzer {
21
22//**********************************************************************
23template<typename EvalT, typename Traits>
26 const Teuchos::ParameterList& p) :
27 residual( p.get<std::string>("Residual Name"),
28 p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
29 scalar( p.get<std::string>("Value Name"),
30 p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar),
31 basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
32{
33 Teuchos::RCP<Teuchos::ParameterList> valid_params = this->getValidParameters();
34 p.validateParameters(*valid_params);
35
36 Teuchos::RCP<const PureBasis> basis
37 = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
38
39 // Verify that this basis supports the gradient operation
40 TEUCHOS_TEST_FOR_EXCEPTION(!basis->isScalarBasis(),std::logic_error,
41 "Integrator_TransientBasisTimesScalar: Basis of type \"" << basis->name() << "\" is not a "
42 "scalar basis");
43
44 this->addEvaluatedField(residual);
45 this->addDependentField(scalar);
46
47 multiplier = p.get<double>("Multiplier");
48
49
50 if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
51 const std::vector<std::string>& field_multiplier_names =
52 *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
53
54 for (std::vector<std::string>::const_iterator name =
55 field_multiplier_names.begin();
56 name != field_multiplier_names.end(); ++name) {
57 PHX::MDField<const ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
58 field_multipliers.push_back(tmp_field);
59 }
60 }
61
62 for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
63 field != field_multipliers.end(); ++field)
64 this->addDependentField(*field);
65
66 std::string n = "Integrator_TransientBasisTimesScalar: " + residual.fieldTag().name();
67 this->setName(n);
68}
69
70//**********************************************************************
71template<typename EvalT, typename Traits>
72void
75 typename Traits::SetupData sd,
77{
78 this->utils.setFieldData(residual,fm);
79 this->utils.setFieldData(scalar,fm);
80
81 num_nodes = residual.extent(1);
82 num_qp = scalar.extent(1);
83
84 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
85
86 tmp = Kokkos::createDynRankView(residual.get_static_view(),"tmp",scalar.extent(0), num_qp);
87}
88
89//**********************************************************************
90template<typename EvalT, typename Traits>
91void
94 typename Traits::EvalData workset)
95{
96 if (workset.evaluate_transient_terms) {
97
98 // for (int i=0; i < residual.size(); ++i)
99 // residual[i] = 0.0;
100
101 Kokkos::deep_copy (residual.get_static_view(), ScalarT(0.0));
102
103 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
104 for (std::size_t qp = 0; qp < num_qp; ++qp) {
105 tmp(cell,qp) = multiplier * scalar(cell,qp);
106 for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
107 field != field_multipliers.end(); ++field)
108 tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
109 }
110 }
111
112 if(workset.num_cells>0)
113 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
114 integrate<ScalarT>(residual.get_view(),
115 tmp,
116 (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
117 }
118}
119
120//**********************************************************************
121template<typename EvalT, typename TRAITS>
122Teuchos::RCP<Teuchos::ParameterList>
124{
125 Teuchos::RCP<Teuchos::ParameterList> p = Teuchos::rcp(new Teuchos::ParameterList);
126 p->set<std::string>("Residual Name", "?");
127 p->set<std::string>("Value Name", "?");
128 Teuchos::RCP<panzer::BasisIRLayout> basis;
129 p->set("Basis", basis);
130 Teuchos::RCP<panzer::IntegrationRule> ir;
131 p->set("IR", ir);
132 p->set<double>("Multiplier", 1.0);
133 Teuchos::RCP<const std::vector<std::string> > fms;
134 p->set("Field Multipliers", fms);
135 return p;
136}
137
138//**********************************************************************
139
140}
141
142#endif
143
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int num_cells
DEPRECATED - use: numCells()
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< const std::vector< panzer::Workset > > worksets_