Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ConstantVector_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_CONSTANT_VECTOR_IMPL_HPP
12#define PANZER_CONSTANT_VECTOR_IMPL_HPP
13
14namespace panzer {
15
16//**********************************************************************
17template<typename EvalT, typename Traits>
20 const Teuchos::ParameterList& p) :
21 vec_(p.get<std::string>("Name"),
22 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout") )
23{
24 this->addEvaluatedField(vec_);
25
26 // Make this unshared so that it is not overwritten
27 this->addUnsharedField(vec_.fieldTag().clone());
28
29 const int dim = vec_.fieldTag().dataLayout().extent(2);
30
31 vals_ = Kokkos::View<double*>("ConstantVector::vals",dim);
32 auto vals_host = Kokkos::create_mirror_view(Kokkos::HostSpace(),vals_);
33
34 vals_host(0) = p.get<double>("Value X");
35 if(dim>1)
36 vals_host(1) = p.get<double>("Value Y");
37 if(dim>2)
38 vals_host(2) = p.get<double>("Value Z");
39
40 Kokkos::deep_copy(vals_,vals_host);
41
42 std::string n = "ConstantVector: " + vec_.fieldTag().name();
43 this->setName(n);
44}
45
46//**********************************************************************
47template<typename EvalT, typename Traits>
48void
50postRegistrationSetup(typename Traits::SetupData /* worksets */,
52{
53 auto vals = this->vals_;
54 auto vec = this->vec_;
55 Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<3>> policy({0,0,0},{static_cast<int64_t>(vec.extent(0)),
56 static_cast<int64_t>(vec.extent(1)),static_cast<int64_t>(vec.extent(2))});
57 Kokkos::parallel_for("panzer::ConstantVector",policy,KOKKOS_LAMBDA(const int c, const int p, const int d){
58 vec(c,p,d) = vals(d);
59 });
60}
61
62//**********************************************************************
63template<typename EvalT, typename Traits>
64void
68
69//**********************************************************************
70
71}
72
73#endif
Kokkos::View< typename Sacado::ScalarType< ScalarT >::type **, Kokkos::LayoutRight, PHX::Device > vals
Kokkos::View< double * > vals_
PHX::MDField< ScalarT > vec_
ConstantVector(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)