Panzer  Version of the Day
Panzer_DOF_BasisToBasis_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_DOF_BASIS_TO_BASIS_IMPL_HPP
44 #define PANZER_DOF_BASIS_TO_BASIS_IMPL_HPP
45 
48 #include "Intrepid2_FunctionSpaceTools.hpp"
49 #include "Phalanx_DataLayout_MDALayout.hpp"
50 #include "Intrepid2_Basis.hpp"
51 #include "Teuchos_Assert.hpp"
52 
53 namespace panzer {
54 
55 //**********************************************************************
56 template <typename EvalT, typename TRAITST>
58 DOF_BasisToBasis(const std::string & fieldName,
59  const PureBasis & sourceBasis,
60  const PureBasis & targetBasis)
61 {
62  TEUCHOS_ASSERT(sourceBasis.numCells() == targetBasis.numCells());
63 
64  // **************
65  // Declare fields
66  // **************
67 
68  dof_source_coeff = PHX::MDField<ScalarT>(fieldName,sourceBasis.functional);
69  dof_target_coeff = PHX::MDField<ScalarT>(fieldName,targetBasis.functional);
70 
71  this->addDependentField(dof_source_coeff);
72  this->addEvaluatedField(dof_target_coeff);
73 
74  // **************
75  // Get coordinate points for reference cell on target basis
76  // **************
78  Teuchos::rcp_dynamic_cast<Intrepid2::DofCoordsInterface<Kokkos::DynRankView<double,PHX::Device> > >(targetBasis.getIntrepid2Basis(),true);
79 
80  Kokkos::DynRankView<double,PHX::Device>intrpCoords =
81  Kokkos::DynRankView<double,PHX::Device>("intrpCoords",targetBasis.cardinality(),targetBasis.dimension());
82 
83  coords_interface->getDofCoords(intrpCoords);
84 
85  // **************
86  // Evaluate source basis values at target basis coordinates
87  // **************
88  Kokkos::DynRankView<double,PHX::Device> basisRef =
89  Kokkos::DynRankView<double,PHX::Device>("basisRef",sourceBasis.cardinality(),targetBasis.cardinality());
90 
91  sourceBasis.getIntrepid2Basis()->getValues(basisRef, intrpCoords, Intrepid2::OPERATOR_VALUE);
92 
93  // **************
94  // Copy the reference basis values for all cells in workset
95  // **************
96  basis = Kokkos::DynRankView<double,PHX::Device>("basis",sourceBasis.numCells(),sourceBasis.cardinality(),targetBasis.cardinality());
97  Intrepid2::FunctionSpaceTools::HGRADtransformVALUE<double>(basis,basisRef);
98 
99  std::string n = "DOF_BasisToBasis: " + dof_target_coeff.fieldTag().name();
100  this->setName(n);
101 }
102 
103 //**********************************************************************
104 template <typename EvalT, typename TRAITST>
107 {
108  this->utils.setFieldData(dof_source_coeff,fm);
109  this->utils.setFieldData(dof_target_coeff,fm);
110 }
111 
112 //**********************************************************************
113 template <typename EvalT, typename TRAITST>
114 void DOF_BasisToBasis<EvalT,TRAITST>::evaluateFields(typename TRAITST::EvalData workset)
115 {
116  // Zero out arrays (intrepid does a sum!)
117  dof_target_coeff.deep_copy(ScalarT(0.0));
118 
119  if(workset.num_cells>0) {
120 
121  // evaluate function at specified points
122  Intrepid2::FunctionSpaceTools::evaluate<ScalarT>(dof_target_coeff,dof_source_coeff,basis);
123  }
124 }
125 
126 //**********************************************************************
127 
128 }
129 
130 #endif
void postRegistrationSetup(typename TRAITST::SetupData d, PHX::FieldManager< TRAITST > &vm)
DOF_BasisToBasis(const std::string &fieldName, const PureBasis &sourceBasis, const PureBasis &targetBasis)
Ctor.
int dimension() const
Returns the dimension of the basis from the topology.
Teuchos::RCP< Intrepid2::Basis< double, Kokkos::DynRankView< double, PHX::Device > > > getIntrepid2Basis() const
int numCells() const
Returns the number of cells in the data layouts.
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
void evaluateFields(typename TRAITST::EvalData workset)
Description and data layouts associated with a particular basis.
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
int cardinality() const
Returns the number of basis coefficients.