Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_DOF_PointField_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_POINT_FIELD_DECL_HPP
44 #define PANZER_DOF_POINT_FIELD_DECL_HPP
45 
48 
49 #include "Intrepid2_FunctionSpaceTools.hpp"
50 
51 #include "Phalanx_DataLayout_MDALayout.hpp"
52 
53 namespace panzer {
54 
55 //**********************************************************************
56 template <typename EvalT, typename TRAITST>
57 void DOF_PointField<EvalT,TRAITST>::initialize(const std::string & fieldName,
58  const PureBasis & fieldBasis,
59  const std::string & coordinateName,
60  const Teuchos::RCP<PHX::DataLayout> & coordLayout,
61  const Teuchos::RCP<PHX::DataLayout> & quadLayout,
62  const std::string & postfixFieldName)
63 {
64  intrepidBasis = fieldBasis.getIntrepid2Basis();
65 
66  int cellCount = fieldBasis.functional->extent(0);
67  int coeffCount = fieldBasis.functional->extent(1);
68  int pointCount = coordLayout->extent(0);
69  int dimCount = coordLayout->extent(1);
70 
71  Teuchos::RCP<PHX::DataLayout> basisLayout = fieldBasis.functional;
72 
73  coordinates = PHX::MDField<const ScalarT,Point,Dim>(coordinateName,coordLayout);
74  dof_coeff = PHX::MDField<const ScalarT>(fieldName,basisLayout);
75  dof_field = PHX::MDField<ScalarT>(fieldName+postfixFieldName,quadLayout);
76 
77  this->addDependentField(coordinates);
78  this->addDependentField(dof_coeff);
79  this->addEvaluatedField(dof_field);
80 
81  // build data storage for temporary conversion
82  basisRef = Kokkos::DynRankView<double,PHX::Device>("basisRef",coeffCount,pointCount);
83  basis = Kokkos::DynRankView<double,PHX::Device>("basis",cellCount,coeffCount,pointCount);
84  intrpCoords = Kokkos::DynRankView<double,PHX::Device>("intrpCoords",pointCount,dimCount);
85 
86  std::string n = "DOF_PointField: " + dof_field.fieldTag().name();
87  this->setName(n);
88 }
89 
90 //**********************************************************************
91 template <typename EvalT, typename TRAITST>
92 void DOF_PointField<EvalT,TRAITST>::evaluateFields(typename TRAITST::EvalData workset)
93 {
94  // Zero out arrays (intrepid does a sum! 1/17/2012)
95  dof_field.deep_copy(ScalarT(0.0));
96 
97  // copy coordinates
98  auto l_intrpCoords = PHX::as_view(intrpCoords);
99  auto l_coordinates = coordinates.get_static_view();
100  Kokkos::parallel_for("DOF PointFields", l_coordinates.extent_int(0), KOKKOS_LAMBDA (int i) {
101  for (int j = 0; j < l_coordinates.extent_int(1); ++j)
102  l_intrpCoords(i,j) = Sacado::scalarValue(l_coordinates(i,j));
103  });
104 
105  if(workset.num_cells>0) {
106  // evaluate at reference points
107  intrepidBasis->getValues(basisRef, intrpCoords, Intrepid2::OPERATOR_VALUE);
108 
109  // transfer reference basis values to physical frame values
110  Intrepid2::FunctionSpaceTools<PHX::exec_space>::
111  HGRADtransformVALUE(basis,basisRef);
112 
113  // evaluate function at specified points
114  Intrepid2::FunctionSpaceTools<PHX::exec_space>::
115  evaluate(dof_field.get_view(),dof_coeff.get_view(),basis);
116  }
117  Kokkos::fence();
118 }
119 
120 //**********************************************************************
121 
122 }
123 
124 #endif
void initialize(const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout, const std::string &postfixFieldName)
Convenience initialization routine, see constructor above.
void evaluateFields(typename TRAITST::EvalData workset)
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
Description and data layouts associated with a particular basis.
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt; or &lt;Cell,Basis&gt;