Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_VectorToScalar_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_VECTOR_TO_SCALAR_IMPL_HPP
12#define PANZER_VECTOR_TO_SCALAR_IMPL_HPP
13
14#include <string>
15
16namespace panzer {
17
18//**********************************************************************
19template<typename EvalT, typename Traits>
22 const Teuchos::ParameterList& p)
23{
24 Teuchos::RCP<PHX::DataLayout> scalar_dl =
25 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Scalar");
26
27 Teuchos::RCP<PHX::DataLayout> vector_dl =
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<ScalarT,Cell,Point>(scalar_names[i], scalar_dl);
37
38 vector_field =
39 PHX::MDField<const ScalarT,Cell,Point,Dim>(p.get<std::string>
40 ("Vector Name"), vector_dl);
41
42 this->addDependentField(vector_field);
43
44 for (std::size_t i=0; i < scalar_fields.size(); ++i)
45 this->addEvaluatedField(scalar_fields[i]);
46
47 std::string n = "VectorToScalar: " + vector_field.fieldTag().name();
48 this->setName(n);
49}
50
51//**********************************************************************
52
53template<typename EvalT, typename Traits> \
54VectorToScalar<EvalT,Traits>::
55VectorToScalar(const PHX::FieldTag & input,
56 const std::vector<PHX::Tag<ScalarT>> & output)
57{
58 // setup the fields
59 vector_field = input;
60
61 scalar_fields.resize(output.size());
62 for(std::size_t i=0;i<output.size();i++)
63 scalar_fields[i] = output[i];
64
65 // add dependent/evaluate fields
66 this->addDependentField(vector_field);
67
68 for (std::size_t i=0; i < scalar_fields.size(); ++i)
69 this->addEvaluatedField(scalar_fields[i]);
70
71 // name array
72 std::string n = "VectorToScalar: " + vector_field.fieldTag().name();
73 this->setName(n);
74}
75
76//**********************************************************************
77template<typename EvalT, typename Traits>
78void
81 typename Traits::EvalData workset)
82{
83
84 // Iteration bounds
85 const int num_points = vector_field.extent_int(1);
86 const int num_scalars = std::min(static_cast<int>(scalar_fields.size()),
87 vector_field.extent_int(2));
88
89 // Need local copies for cuda, *this is not usable
90 auto local_vector_field = vector_field;
91
92 // We parallelize over each scalar field
93 for (int sc = 0; sc < num_scalars; ++sc) {
94 auto local_scalar_field = scalar_fields[sc];
95 Kokkos::parallel_for(workset.num_cells, KOKKOS_LAMBDA (const int cell) {
96 for (int pt = 0; pt < num_points; ++pt)
97 local_scalar_field(cell,pt) = local_vector_field(cell,pt,sc);
98 });
99 }
100
101 // If there are remaining fields, set them to zero
102 for(unsigned int sc = num_scalars; sc < scalar_fields.size(); ++sc)
103 Kokkos::deep_copy(scalar_fields[sc].get_static_view(), 0.);
104}
105
106//**********************************************************************
107
108}
109
110#endif
void evaluateFields(typename Traits::EvalData d)
VectorToScalar(const Teuchos::ParameterList &p)
int num_cells
DEPRECATED - use: numCells()