Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherBasisCoordinates_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_GATHER_BASIS_COORDINATES_IMPL_HPP
12 #define PANZER_GATHER_BASIS_COORDINATES_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 
17 #include "Panzer_GlobalIndexer.hpp"
18 #include "Panzer_PureBasis.hpp"
20 
21 #include "Teuchos_FancyOStream.hpp"
22 
23 template<typename EvalT,typename TRAITS>
24 std::string
26 fieldName(const std::string & basisName)
27 {
28  std::stringstream ss;
29  ss << "Basis_" << basisName << " BasisCoordinates";
30  return ss.str();
31 }
32 
33 template<typename EvalT,typename TRAITS>
36 {
37  basisName_ = basis.name();
38 
39  basisCoordinates_ = PHX::MDField<ScalarT,Cell,BASIS,Dim>(fieldName(basisName_),basis.coordinates);
40 
41  this->addEvaluatedField(basisCoordinates_);
42  this->addUnsharedField(Teuchos::rcp_const_cast<PHX::FieldTag>(basisCoordinates_.fieldTagPtr()));
43 
44  this->setName("GatherBasisCoordinates: "+fieldName(basisName_));
45 }
46 
47 // **********************************************************************
48 template<typename EvalT,typename TRAITS>
50 postRegistrationSetup(typename TRAITS::SetupData sd,
51  PHX::FieldManager<TRAITS>& /* fm */)
52 {
53  basisIndex_ = panzer::getPureBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
54  // make sure to zero out derivative array as we only set the value for AD types in the loop below
55  Kokkos::deep_copy(basisCoordinates_.get_static_view(),0.0);
56 }
57 
58 // **********************************************************************
59 template<typename EvalT,typename TRAITS>
61 evaluateFields(typename TRAITS::EvalData workset)
62 {
63  // const Kokkos::DynRankView<double,PHX::Device> & basisCoords = this->wda(workset).bases[basisIndex_]->basis_coordinates;
64  const Teuchos::RCP<const BasisValues2<double> > bv = this->wda(workset).bases[basisIndex_];
65 
66  // just copy the array
67  auto d_basisCoordinates = basisCoordinates_.get_static_view();
68  auto s_basis_coordinates = bv->basis_coordinates.get_static_view();
69 
70  Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<3>> policy({0,0,0},{int(workset.num_cells),s_basis_coordinates.extent_int(1),s_basis_coordinates.extent_int(2)});
71  Kokkos::parallel_for("GatherBasisCoords",policy, KOKKOS_LAMBDA(const int i, const int j, const int k) {
72  auto d_basisCoordinates_tmp = d_basisCoordinates;
73  auto s_basis_coordinates_tmp = s_basis_coordinates;
74  if constexpr(Sacado::IsADType<typename EvalT::ScalarT>::value) {
75  d_basisCoordinates_tmp(i,j,k).val() = s_basis_coordinates_tmp(i,j,k);
76  }
77  else {
78  d_basisCoordinates_tmp(i,j,k) = s_basis_coordinates_tmp(i,j,k);
79  }
80  });
81  Kokkos::fence();
82 }
83 
84 #endif
std::string name() const
A unique key that is the combination of the basis type and basis order.
void evaluateFields(typename TRAITS::EvalData d)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)
Teuchos::RCP< PHX::DataLayout > coordinates
&lt;Cell,Basis,Dim&gt;
std::vector< std::string >::size_type getPureBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular PureBasis name.
Description and data layouts associated with a particular basis.
static std::string fieldName(const std::string &basisName)