Panzer Version of the Day
Loading...
Searching...
No Matches
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
16namespace panzer {
17
18//**********************************************************************
19template<typename EvalT, typename Traits>
21ScalarToVector(const Teuchos::ParameterList& p)
22{
23 auto scalar_dl = p.get<Teuchos::RCP<PHX::DataLayout>>("Data Layout Scalar");
24 auto vector_dl = p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
25 auto scalar_names = p.get<Teuchos::RCP<const std::vector<std::string>>>("Scalar Names");
26 auto vector_name = p.get<std::string>("Vector Name");
27
28 scalar_fields_.resize(scalar_names->size());
29 for (size_t i=0; i < scalar_names->size(); ++i)
30 scalar_fields_[i] = PHX::MDField<const ScalarT,Cell,Point>((*scalar_names)[i], scalar_dl);
31
32 vector_field_ = PHX::MDField<ScalarT,Cell,Point,Dim>(vector_name, vector_dl);
33
34 for (std::size_t i=0; i < scalar_fields_.size(); ++i)
35 this->addDependentField(scalar_fields_[i]);
36
37 this->addEvaluatedField(vector_field_);
38
39 std::string n = "ScalarToVector: " + vector_field_.fieldTag().name();
40 this->setName(n);
41}
42
43//**********************************************************************
44
45template<typename EvalT, typename Traits> \
46ScalarToVector<EvalT,Traits>::
47ScalarToVector(const std::vector<PHX::Tag<ScalarT>> & input,
48 const PHX::FieldTag & output)
49{
50 // setup the fields
51 vector_field_ = output;
52
53 scalar_fields_.resize(input.size());
54 for(std::size_t i=0;i<input.size();i++)
55 scalar_fields_[i] = input[i];
56
57 // add dependent/evaluate fields
58 this->addEvaluatedField(vector_field_);
59
60 for (std::size_t i=0; i < scalar_fields_.size(); ++i)
61 this->addDependentField(scalar_fields_[i]);
62
63 // name array
64 std::string n = "ScalarToVector: " + vector_field_.fieldTag().name();
65 this->setName(n);
66}
67
68//**********************************************************************
69template<typename EvalT, typename Traits>
70void
73 typename Traits::SetupData /* worksets */,
75{
76 scalar_fields_vov_.initialize("Scalar Fields VoV",scalar_fields_.size());
77 for (size_t i=0; i < scalar_fields_.size(); ++i)
78 scalar_fields_vov_.addView(scalar_fields_[i].get_static_view(),i);
79
80 scalar_fields_vov_.syncHostToDevice();
81}
82
83//**********************************************************************
84template<typename EvalT, typename Traits>
85void
88 typename Traits::EvalData workset)
89{
90 const int num_points = vector_field_.extent_int(1);
91 const int num_vector_scalars = vector_field_.extent_int(2);
92 auto vec = vector_field_;
93 auto scalars = scalar_fields_vov_.getViewDevice();
94
95 // Corner case: user can supply fewer scalar fields than the
96 // dimension on vector. If so, fill missing fields with zeroes.
97 const int num_scalars = std::min(static_cast<int>(scalars.extent(0)),num_vector_scalars);
98
99 Kokkos::parallel_for (workset.num_cells,KOKKOS_LAMBDA (const int cell){
100 for (int pt = 0; pt < num_points; ++pt) {
101 for (int sc = 0; sc < num_scalars; ++sc) {
102 vec(cell,pt,sc) = scalars(sc)(cell,pt);
103 }
104
105 // Missing scalars are filled with zero
106 for(int sc = num_scalars; sc < num_vector_scalars; ++sc) {
107 vec(cell,pt,sc) = ScalarT(0.0);
108 }
109 }
110 });
111}
112
113//**********************************************************************
114
115}
116
117#endif
ScalarToVector(const Teuchos::ParameterList &p)
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int num_cells
DEPRECATED - use: numCells()