Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ScalarToVector_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_SCALAR_TO_VECTOR_IMPL_HPP
12 #define PANZER_SCALAR_TO_VECTOR_IMPL_HPP
13 
14 #include <string>
15 
16 namespace panzer {
17 
18 //**********************************************************************
19 template<typename EvalT, typename Traits>
22  const Teuchos::ParameterList& p)
23 {
25  p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Scalar");
26 
28  p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
29 
30  const std::vector<std::string>& scalar_names =
31  *(p.get< Teuchos::RCP<const std::vector<std::string> > >("Scalar Names"));
32 
33  scalar_fields.resize(scalar_names.size());
34  for (std::size_t i=0; i < scalar_names.size(); ++i)
35  scalar_fields[i] =
36  PHX::MDField<const ScalarT,Cell,Point>(scalar_names[i], scalar_dl);
37 
38  vector_field =
40  ("Vector Name"), vector_dl);
41 
42  for (std::size_t i=0; i < scalar_fields.size(); ++i)
43  this->addDependentField(scalar_fields[i]);
44 
45  this->addEvaluatedField(vector_field);
46 
47  std::string n = "ScalarToVector: " + vector_field.fieldTag().name();
48  this->setName(n);
49 }
50 
51 //**********************************************************************
52 
53 template<typename EvalT, typename Traits> \
54 ScalarToVector<EvalT,Traits>::
55 ScalarToVector(const std::vector<PHX::Tag<ScalarT>> & input,
56  const PHX::FieldTag & output)
57 {
58  // setup the fields
59  vector_field = output;
60 
61  scalar_fields.resize(input.size());
62  for(std::size_t i=0;i<input.size();i++)
63  scalar_fields[i] = input[i];
64 
65  // add dependent/evaluate fields
66  this->addEvaluatedField(vector_field);
67 
68  for (std::size_t i=0; i < scalar_fields.size(); ++i)
69  this->addDependentField(scalar_fields[i]);
70 
71  // name array
72  std::string n = "ScalarToVector: " + vector_field.fieldTag().name();
73  this->setName(n);
74 }
75 
76 //**********************************************************************
77 template<typename EvalT, typename Traits>
78 void
81  typename Traits::SetupData /* worksets */,
83 {
84  // Convert std::vector to PHX::View for use on device
85  internal_scalar_fields = PHX::View<KokkosScalarFields_t*>("ScalarToVector::internal_scalar_fields", scalar_fields.size());
86  for (std::size_t i=0; i < scalar_fields.size(); ++i)
87  internal_scalar_fields(i) = scalar_fields[i].get_static_view();
88 }
89 
90 //**********************************************************************
91 template<typename EvalT, typename Traits>
92 void
95  typename Traits::EvalData workset)
96 {
97 
98  using Scalar = typename EvalT::ScalarT;
99 
100  // Iteration bounds
101  const int num_points = vector_field.extent_int(1);
102  const int num_vector_scalars = vector_field.extent_int(2);
103  const int num_scalars = std::min(internal_scalar_fields.extent_int(0),
104  num_vector_scalars);
105 
106  // Local copies to prevent passing (*this) to device code
107  auto vector = vector_field;
108  auto scalars = internal_scalar_fields;
109 
110  // Loop over cells, points, scalars
111  Kokkos::parallel_for (workset.num_cells,KOKKOS_LAMBDA (const int cell){
112  for (int pt = 0; pt < num_points; ++pt){
113 
114  // Copy over scalars
115  for (int sc = 0; sc < num_scalars; ++sc)
116  vector(cell,pt,sc) = scalars(sc)(cell,pt);
117 
118  // Missing scalars are filled with zero
119  for(int sc = num_scalars; sc < num_vector_scalars; ++sc)
120  vector(cell,pt,sc) = Scalar(0);
121  }
122  });
123 
124 }
125 
126 //**********************************************************************
127 
128 }
129 
130 #endif
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
ScalarToVector(const Teuchos::ParameterList &p)