Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Parameter_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_PARAMETER_IMPL_HPP
12 #define PANZER_PARAMETER_IMPL_HPP
13 
14 #include "PanzerDiscFE_config.hpp"
17 #include <cstddef>
18 #include <string>
19 #include <vector>
20 
21 namespace panzer {
22 
23 //**********************************************************************
24 template<typename EvalT, typename TRAITS>
26 Parameter(const std::string parameter_name,
27  const std::string field_name,
28  const Teuchos::RCP<PHX::DataLayout>& data_layout,
29  panzer::ParamLib& param_lib)
30 {
31  target_field = PHX::MDField<ScalarT, Cell, Point>(field_name, data_layout);
32 
33  this->addEvaluatedField(target_field);
34 
35  //param = panzer::accessScalarParameter<EvalT>(parameter_name,param_lib);
36  param = panzer::createAndRegisterScalarParameter<EvalT>(parameter_name,param_lib);
37  // no initialization, this will be done by someone else (possibly the ME) later
38 
39  std::string n = "Parameter Evaluator";
40  this->setName(n);
41 }
42 
43 //**********************************************************************
44 template<typename EvalT, typename TRAITS>
46 evaluateFields(typename TRAITS::EvalData workset)
47 {
48  //std::cout << "Parameter::evalauteFields() ParamValue = " << param->getValue() << std::endl;
49  auto param_val = param->getValue();
50  auto target_field_v = target_field.get_static_view();
51  auto target_field_h = Kokkos::create_mirror_view(target_field_v);
52 
53  for (int cell=0; cell < workset.num_cells; ++cell) {
54  for (std::size_t pt=0; pt<target_field_v.extent(1); ++pt)
55  target_field_h(cell,pt) = param_val;
56  }
57  Kokkos::deep_copy(target_field_v, target_field_h);
58 
59 }
60 
61 //**********************************************************************
62 
63 }
64 
65 #endif
Sacado::ScalarParameterLibrary< panzer::EvaluationTraits > ParamLib
Parameter(const std::string parameter_name, const std::string field_name, const Teuchos::RCP< PHX::DataLayout > &data_layout, panzer::ParamLib &param_lib)
void evaluateFields(typename TRAITS::EvalData ud)