Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Dirichlet_Residual_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_DIRICHLET_RESIDUAL_IMPL_HPP
12#define PANZER_DIRICHLET_RESIDUAL_IMPL_HPP
13
14#include <cstddef>
15#include <string>
16#include <vector>
17
18namespace panzer {
19
20//**********************************************************************
21template<typename EvalT, typename Traits>
24 const Teuchos::ParameterList& p)
25{
26 std::string residual_name = p.get<std::string>("Residual Name");
27 std::string dof_name = p.get<std::string>("DOF Name");
28 std::string value_name = p.get<std::string>("Value Name");
29
30 Teuchos::RCP<PHX::DataLayout> data_layout =
31 p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout");
32
33 residual = PHX::MDField<ScalarT>(residual_name, data_layout);
34 dof = PHX::MDField<const ScalarT>(dof_name, data_layout);
35 value = PHX::MDField<const ScalarT>(value_name, data_layout);
36
37 this->addEvaluatedField(residual);
38 this->addDependentField(dof);
39 this->addDependentField(value);
40
41 std::string n = "Dirichlet Residual Evaluator";
42 this->setName(n);
43}
44
45//**********************************************************************
46template<typename EvalT, typename Traits>
47void
50 typename Traits::SetupData /* worksets */,
52{
53 cell_data_size = residual.fieldTag().dataLayout().extent(1);
54}
55
56//**********************************************************************
57template<typename EvalT, typename Traits>
58void
61 typename Traits::EvalData workset)
62{
63 auto residual_v = residual.get_static_view();
64 auto dof_v = dof.get_static_view();
65 auto value_v = value.get_static_view();
66 auto local_cell_data_size = cell_data_size;
67 Kokkos::parallel_for (workset.num_cells, KOKKOS_LAMBDA (index_t i) {
68 for (std::size_t j = 0; j < local_cell_data_size; ++j)
69 residual_v(i,j)=dof_v(i,j)-value_v(i,j);
70 });
71}
72
73//**********************************************************************
74
75}
76
77#endif
void evaluateFields(typename Traits::EvalData d)
DirichletResidual(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int num_cells
DEPRECATED - use: numCells()