Panzer  Version of the Day
Panzer_ScalarToVector_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_SCALAR_TO_VECTOR_IMPL_HPP
44 #define PANZER_SCALAR_TO_VECTOR_IMPL_HPP
45 
46 #include <string>
47 
48 namespace panzer {
49 
50 //**********************************************************************
51 PHX_EVALUATOR_CTOR(ScalarToVector,p)
52 {
53  Teuchos::RCP<PHX::DataLayout> scalar_dl =
54  p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Scalar");
55 
56  Teuchos::RCP<PHX::DataLayout> vector_dl =
57  p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
58 
59  const std::vector<std::string>& scalar_names =
60  *(p.get< Teuchos::RCP<const std::vector<std::string> > >("Scalar Names"));
61 
62  scalar_fields.resize(scalar_names.size());
63  for (std::size_t i=0; i < scalar_names.size(); ++i)
64  scalar_fields[i] =
65  PHX::MDField<ScalarT,Cell,Point>(scalar_names[i], scalar_dl);
66 
67  vector_field =
68  PHX::MDField<ScalarT,Cell,Point,Dim>(p.get<std::string>
69  ("Vector Name"), vector_dl);
70 
71  for (std::size_t i=0; i < scalar_fields.size(); ++i)
72  this->addDependentField(scalar_fields[i]);
73 
74  this->addEvaluatedField(vector_field);
75 
76  std::string n = "ScalarToVector: " + vector_field.fieldTag().name();
77  this->setName(n);
78 }
79 
80 //**********************************************************************
81 PHX_POST_REGISTRATION_SETUP(ScalarToVector,worksets,fm)
82 {
83  internal_scalar_fields = Kokkos::View<KokkosScalarFields_t*>("ScalarToVector::internal_scalar_fields", scalar_fields.size());
84  for (std::size_t i=0; i < scalar_fields.size(); ++i) {
85  this->utils.setFieldData(scalar_fields[i],fm);
86  internal_scalar_fields(i) = scalar_fields[i].get_static_view();
87  }
88 
89  this->utils.setFieldData(vector_field,fm);
90 }
91 
92 //**********************************************************************
93 
94 template<typename EvalT, typename TRAITS>
95 KOKKOS_INLINE_FUNCTION
96 void ScalarToVector<EvalT, TRAITS>::operator()(const size_t &cell) const {
97  typedef typename PHX::MDField<ScalarT,Cell,Point>::size_type size_type;
98  // Loop over points
99  for (size_type pt = 0; pt < vector_field.dimension(1); ++pt) {
100  // Loop over scalars
101  for (std::size_t sc = 0; sc < internal_scalar_fields.dimension(0); ++sc) {
102  vector_field(cell,pt,sc) = internal_scalar_fields(sc)(cell,pt);
103 
104  }
105  }
106 
107 }
108 //**********************************************************************
109 PHX_EVALUATE_FIELDS(ScalarToVector,workset)
110 {
111 
112  // Loop over cells
113  Kokkos::parallel_for (workset.num_cells, (*this));
114 }
115 
116 //**********************************************************************
117 
118 }
119 
120 #endif
std::vector< PHX::MDField< ScalarT, Cell, Point > > scalar_fields
Interpolates basis DOF values to IP DOF values.
void operator()(const FieldMultTag< NUM_FIELD_MULT > &, const std::size_t &cell) const
PHX::MDField< ScalarT, Cell, Point, Dim > vector_field
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Kokkos::View< KokkosScalarFields_t * > internal_scalar_fields
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)