43 #ifndef PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP 44 #define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP 46 #include "Intrepid2_FunctionSpaceTools.hpp" 50 #include "Kokkos_ViewFactory.hpp" 62 Teuchos::RCP<const BasisIRLayout> basis_layout = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis").getConst();
63 Teuchos::RCP<const panzer::IntegrationRule> ir = p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR").getConst();
64 Teuchos::RCP<PHX::DataLayout> vector_dl = p.get< Teuchos::RCP<PHX::DataLayout> >(
"Data Layout Vector");
66 Teuchos::RCP<const PureBasis>
basis = basis_layout->getBasis();
70 TEUCHOS_TEST_FOR_EXCEPTION(!
basis->supportsGrad(),std::logic_error,
"Integrator_GradBasisCrossVector: Basis of type \"" <<
basis->name() <<
"\" does not support GRAD");
74 if (p.isType<Teuchos::RCP<
const std::vector<std::string> > >(
"Residual Names")){
75 const std::vector<std::string> & names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >(
"Residual Names"));
76 for(
const auto & name : names){
77 PHX::MDField<ScalarT,Cell,BASIS> res(name, basis_layout->functional);
84 TEUCHOS_TEST_FOR_EXCEPTION(
_num_dims != vector_dl->dimension(2),std::logic_error,
"Vector must be same length as number of residuals.");
87 _vector = PHX::MDField<const ScalarT,Cell,IP,Dim>(p.get<std::string>(
"Vector Name"), vector_dl);
93 TEUCHOS_TEST_FOR_EXCEPTION(
_num_grad_dims >
_num_dims,std::logic_error,
"Vector size must have at least as many components as there are dimensions in the mesh.");
99 if (p.isType<Teuchos::RCP<
const std::vector<std::string> > >(
"Field Multipliers")){
100 const std::vector<std::string> & field_multiplier_names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >(
"Field Multipliers"));
101 for (
const std::string & name : field_multiplier_names){
102 _field_multipliers.push_back(PHX::MDField<const ScalarT,Cell,IP>(name, ir->dl_scalar));
112 this->addDependentField(
_vector);
114 this->addDependentField(
field);
118 this->setName(
"Integrator_GradBasisCrossVector: " +
_vector.fieldTag().name());
126 this->utils.setFieldData(
residual,fm);
129 this->utils.setFieldData(
_vector,fm);
132 this->utils.setFieldData(
field,fm);
150 Kokkos::deep_copy(
_residuals[i].get_static_view(), ScalarT(0.0));
159 for (
int i=0; i <
_vector.extent_int(0); ++i)
160 for (
int j=0; j <
_vector.extent_int(1); ++j)
161 for (
int k=0; k <
_vector.extent_int(2); ++k)
166 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
168 const ScalarT & tmpVar = field_data(cell,qp);
169 for (std::size_t dim = 0; dim <
_num_dims; ++dim){
170 _tmp(cell,qp,dim) *= tmpVar;
187 for (index_t cell = 0; cell < workset.num_cells; ++cell){
204 for (index_t cell = 0; cell < workset.num_cells; ++cell){
210 r_2 +=
_tmp(cell,qp,0)*bv.
weighted_grad_basis(cell,
basis,qp,1) -
_tmp(cell,qp,1)*bv.
weighted_grad_basis(cell,
basis,qp,0);
222 for (index_t cell = 0; cell < workset.num_cells; ++cell){
226 r_0 +=
_tmp(cell,qp,1)*bv.
weighted_grad_basis(cell,
basis,qp,2) -
_tmp(cell,qp,2)*bv.
weighted_grad_basis(cell,
basis,qp,1);
227 r_1 +=
_tmp(cell,qp,2)*bv.
weighted_grad_basis(cell,
basis,qp,0) -
_tmp(cell,qp,0)*bv.
weighted_grad_basis(cell,
basis,qp,2);
228 r_2 +=
_tmp(cell,qp,0)*bv.
weighted_grad_basis(cell,
basis,qp,1) -
_tmp(cell,qp,1)*bv.
weighted_grad_basis(cell,
basis,qp,0);
PHX::MDField< ScalarT > residual
Evaluates a Dirichlet BC residual corresponding to a field value.
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::size_t _num_grad_dims
Kokkos::DynRankView< ScalarT, PHX::Device > _tmp
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > _residuals
PHX::MDField< const ScalarT, Cell, IP, Dim > _vector
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.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
Array_CellBasisIPDim weighted_grad_basis
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
std::size_t _num_basis_nodes
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)
std::size_t _num_quadrature_points
std::vector< PHX::MDField< const ScalarT, Cell, IP > > _field_multipliers