Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DOF_PointField_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_DOF_POINT_FIELD_DECL_HPP
12#define PANZER_DOF_POINT_FIELD_DECL_HPP
13
16
17#include "Intrepid2_FunctionSpaceTools.hpp"
18
19#include "Phalanx_DataLayout_MDALayout.hpp"
20
21namespace panzer {
22
23//**********************************************************************
24template <typename EvalT, typename TRAITST>
25void DOF_PointField<EvalT,TRAITST>::initialize(const std::string & fieldName,
26 const PureBasis & fieldBasis,
27 const std::string & coordinateName,
28 const Teuchos::RCP<PHX::DataLayout> & coordLayout,
29 const Teuchos::RCP<PHX::DataLayout> & quadLayout,
30 const std::string & postfixFieldName)
31{
32 intrepidBasis = fieldBasis.getIntrepid2Basis();
33
34 int cellCount = fieldBasis.functional->extent(0);
35 int coeffCount = fieldBasis.functional->extent(1);
36 int pointCount = coordLayout->extent(0);
37 int dimCount = coordLayout->extent(1);
38
39 Teuchos::RCP<PHX::DataLayout> basisLayout = fieldBasis.functional;
40
41 coordinates = PHX::MDField<const ScalarT,Point,Dim>(coordinateName,coordLayout);
42 dof_coeff = PHX::MDField<const ScalarT>(fieldName,basisLayout);
43 dof_field = PHX::MDField<ScalarT>(fieldName+postfixFieldName,quadLayout);
44
45 this->addDependentField(coordinates);
46 this->addDependentField(dof_coeff);
47 this->addEvaluatedField(dof_field);
48
49 // build data storage for temporary conversion
50 basisRef = Kokkos::DynRankView<double,PHX::Device>("basisRef",coeffCount,pointCount);
51 basis = Kokkos::DynRankView<double,PHX::Device>("basis",cellCount,coeffCount,pointCount);
52 intrpCoords = Kokkos::DynRankView<double,PHX::Device>("intrpCoords",pointCount,dimCount);
53
54 std::string n = "DOF_PointField: " + dof_field.fieldTag().name();
55 this->setName(n);
56}
57
58//**********************************************************************
59template <typename EvalT, typename TRAITST>
60void DOF_PointField<EvalT,TRAITST>::evaluateFields(typename TRAITST::EvalData workset)
61{
62 // Zero out arrays (intrepid does a sum! 1/17/2012)
63 dof_field.deep_copy(ScalarT(0.0));
64
65 // copy coordinates
66 auto l_intrpCoords = PHX::as_view(intrpCoords);
67 auto l_coordinates = coordinates.get_static_view();
68 Kokkos::parallel_for("DOF PointFields", l_coordinates.extent_int(0), KOKKOS_LAMBDA (int i) {
69 for (int j = 0; j < l_coordinates.extent_int(1); ++j)
70 l_intrpCoords(i,j) = Sacado::scalarValue(l_coordinates(i,j));
71 });
72
73 if(workset.num_cells>0) {
74 // evaluate at reference points
75 intrepidBasis->getValues(basisRef, intrpCoords, Intrepid2::OPERATOR_VALUE);
76
77 // transfer reference basis values to physical frame values
78 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
79 HGRADtransformVALUE(basis,basisRef);
80
81 // evaluate function at specified points
82 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
83 evaluate(dof_field.get_view(),dof_coeff.get_view(),basis);
84 }
85 Kokkos::fence();
86}
87
88//**********************************************************************
89
90}
91
92#endif
void initialize(const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout, const std::string &postfixFieldName)
Convenience initialization routine, see constructor above.
void evaluateFields(typename TRAITST::EvalData workset)
Description and data layouts associated with a particular basis.
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>