Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DOF_BasisToBasis_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_DOF_BASIS_TO_BASIS_IMPL_HPP
12#define PANZER_DOF_BASIS_TO_BASIS_IMPL_HPP
13
16#include "Intrepid2_FunctionSpaceTools.hpp"
17#include "Phalanx_DataLayout_MDALayout.hpp"
18#include "Intrepid2_Basis.hpp"
19#include "Teuchos_Assert.hpp"
20
21namespace panzer {
22
23//**********************************************************************
24template <typename EvalT, typename TRAITST>
26DOF_BasisToBasis(const std::string & fieldName,
27 const PureBasis & sourceBasis,
28 const PureBasis & targetBasis)
29{
30 TEUCHOS_ASSERT(sourceBasis.numCells() == targetBasis.numCells());
31
32 // **************
33 // Declare fields
34 // **************
35 dof_source_coeff = PHX::MDField<const ScalarT>(fieldName,sourceBasis.functional);
36 dof_target_coeff = PHX::MDField<ScalarT>(fieldName,targetBasis.functional);
37
38 this->addDependentField(dof_source_coeff);
39 this->addEvaluatedField(dof_target_coeff);
40
41 // **************
42 // Get coordinate points for reference cell on target basis
43 // **************
44 Kokkos::DynRankView<double,PHX::Device>intrpCoords =
45 Kokkos::DynRankView<double,PHX::Device>("intrpCoords",targetBasis.cardinality(),targetBasis.dimension());
46
47 targetBasis.getIntrepid2Basis<PHX::exec_space,double,double>()->getDofCoords(intrpCoords);
48
49 // **************
50 // Evaluate source basis values at target basis coordinates
51 // **************
52 Kokkos::DynRankView<double,PHX::Device> basisRef =
53 Kokkos::DynRankView<double,PHX::Device>("basisRef",sourceBasis.cardinality(),targetBasis.cardinality());
54
55 sourceBasis.getIntrepid2Basis()->getValues(basisRef, intrpCoords, Intrepid2::OPERATOR_VALUE);
56
57 // **************
58 // Copy the reference basis values for all cells in workset
59 // **************
60 basis = Kokkos::DynRankView<double,PHX::Device>("basis",sourceBasis.numCells(),sourceBasis.cardinality(),targetBasis.cardinality());
61 Intrepid2::FunctionSpaceTools<PHX::exec_space>::HGRADtransformVALUE(basis,basisRef);
62
63 std::string n = "DOF_BasisToBasis: " + dof_target_coeff.fieldTag().name();
64 this->setName(n);
65}
66
67//**********************************************************************
68template <typename EvalT, typename TRAITST>
69void DOF_BasisToBasis<EvalT,TRAITST>::evaluateFields(typename TRAITST::EvalData workset)
70{
71 // Zero out arrays (intrepid does a sum!)
72 dof_target_coeff.deep_copy(ScalarT(0.0));
73
74 if(workset.num_cells>0) {
75
76 // evaluate function at specified points
77 Intrepid2::FunctionSpaceTools<PHX::exec_space>::evaluate(dof_target_coeff.get_view(),
78 dof_source_coeff.get_view(),
79 basis);
80 }
81}
82
83//**********************************************************************
84
85}
86
87#endif
void evaluateFields(typename TRAITST::EvalData workset)
DOF_BasisToBasis(const std::string &fieldName, const PureBasis &sourceBasis, const PureBasis &targetBasis)
Ctor.
Description and data layouts associated with a particular basis.
int numCells() const
Returns the number of cells in the data layouts.
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
int cardinality() const
Returns the number of basis coefficients.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
int dimension() const
Returns the dimension of the basis from the topology.