Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_VectorToScalar_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_VECTOR_TO_SCALAR_IMPL_HPP
44#define PANZER_VECTOR_TO_SCALAR_IMPL_HPP
45
46#include <string>
47
48namespace panzer {
49
50//**********************************************************************
51template<typename EvalT, typename Traits>
54 const Teuchos::ParameterList& p)
55{
56 Teuchos::RCP<PHX::DataLayout> scalar_dl =
57 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Scalar");
58
59 Teuchos::RCP<PHX::DataLayout> vector_dl =
60 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
61
62 const std::vector<std::string>& scalar_names =
63 *(p.get< Teuchos::RCP<const std::vector<std::string> > >("Scalar Names"));
64
65 scalar_fields.resize(scalar_names.size());
66 for (std::size_t i=0; i < scalar_names.size(); ++i)
67 scalar_fields[i] =
68 PHX::MDField<ScalarT,Cell,Point>(scalar_names[i], scalar_dl);
69
70 vector_field =
71 PHX::MDField<const ScalarT,Cell,Point,Dim>(p.get<std::string>
72 ("Vector Name"), vector_dl);
73
74 this->addDependentField(vector_field);
75
76 for (std::size_t i=0; i < scalar_fields.size(); ++i)
77 this->addEvaluatedField(scalar_fields[i]);
78
79 std::string n = "VectorToScalar: " + vector_field.fieldTag().name();
80 this->setName(n);
81}
82
83//**********************************************************************
84
85template<typename EvalT, typename Traits> \
86VectorToScalar<EvalT,Traits>::
87VectorToScalar(const PHX::FieldTag & input,
88 const std::vector<PHX::Tag<ScalarT>> & output)
89{
90 // setup the fields
91 vector_field = input;
92
93 scalar_fields.resize(output.size());
94 for(std::size_t i=0;i<output.size();i++)
95 scalar_fields[i] = output[i];
96
97 // add dependent/evaluate fields
98 this->addDependentField(vector_field);
99
100 for (std::size_t i=0; i < scalar_fields.size(); ++i)
101 this->addEvaluatedField(scalar_fields[i]);
102
103 // name array
104 std::string n = "VectorToScalar: " + vector_field.fieldTag().name();
105 this->setName(n);
106}
107
108//**********************************************************************
109template<typename EvalT, typename Traits>
110void
113 typename Traits::EvalData workset)
114{
115
116 // Iteration bounds
117 const int num_points = vector_field.extent_int(1);
118 const int num_scalars = std::min(static_cast<int>(scalar_fields.size()),
119 vector_field.extent_int(2));
120
121 // Need local copies for cuda, *this is not usable
122 auto local_vector_field = vector_field;
123
124 // We parallelize over each scalar field
125 for (int sc = 0; sc < num_scalars; ++sc) {
126 auto local_scalar_field = scalar_fields[sc];
127 Kokkos::parallel_for(workset.num_cells, KOKKOS_LAMBDA (const int cell) {
128 for (int pt = 0; pt < num_points; ++pt)
129 local_scalar_field(cell,pt) = local_vector_field(cell,pt,sc);
130 });
131 }
132
133 // If there are remaining fields, set them to zero
134 for(unsigned int sc = num_scalars; sc < scalar_fields.size(); ++sc)
135 Kokkos::deep_copy(scalar_fields[sc].get_static_view(), 0.);
136}
137
138//**********************************************************************
139
140}
141
142#endif
void evaluateFields(typename Traits::EvalData d)
VectorToScalar(const Teuchos::ParameterList &p)
int num_cells
DEPRECATED - use: numCells()