Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherIntegrationCoordinates_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_INTEGRATION_COORDINATES_IMPL_HPP
12 #define PANZER_GATHER_INTEGRATION_COORDINATES_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 
19 
20 #include "Teuchos_FancyOStream.hpp"
21 
22 template<typename EvalT,typename TRAITS>
23 std::string
25 fieldName(int degree)
26 {
27  std::stringstream ss;
28  ss << "IR_" << degree << " IntegrationCoordinates";
29  return ss.str();
30 }
31 
32 template<typename EvalT,typename TRAITS>
35 {
36  quadDegree_ = quad.cubature_degree;
37 
38  quadCoordinates_ = PHX::MDField<ScalarT,Cell,Point,Dim>(fieldName(quadDegree_),quad.dl_vector);
39 
40  this->addEvaluatedField(quadCoordinates_);
41  this->addUnsharedField(Teuchos::rcp_const_cast<PHX::FieldTag>(quadCoordinates_.fieldTagPtr()));
42 
43  this->setName("GatherIntegrationCoordinates: "+fieldName(quadDegree_));
44 }
45 
46 // **********************************************************************
47 template<typename EvalT,typename TRAITS>
49 postRegistrationSetup(typename TRAITS::SetupData sd,
50  PHX::FieldManager<TRAITS>& /* fm */)
51 {
52  quadIndex_ = panzer::getIntegrationRuleIndex(quadDegree_, (*sd.worksets_)[0], this->wda);
53  // make sure to zero out derivative array as we only set the value for AD types in the loop below
54  Kokkos::deep_copy(quadCoordinates_.get_static_view(),0.0);
55 }
56 
57 // **********************************************************************
58 template<typename EvalT,typename TRAITS>
60 evaluateFields(typename TRAITS::EvalData workset)
61 {
62  // const Kokkos::DynRankView<double,PHX::Device> & quadCoords = this->wda(workset).int_rules[quadIndex_]->ip_coordinates;
63  const IntegrationValues2<double> & iv = *this->wda(workset).int_rules[quadIndex_];
64  auto s_ip_coordinates = iv.ip_coordinates.get_static_view();
65  auto d_quadCoordinates = quadCoordinates_.get_static_view();
66 
67  // just copy the array
68  Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<3>> policy({0,0,0},{int(workset.num_cells),s_ip_coordinates.extent_int(1),s_ip_coordinates.extent_int(2)});
69  Kokkos::parallel_for("GatherIntegrationCoords", policy, KOKKOS_LAMBDA (const int i, const int j, const int k) {
70  auto s_ip_coordinates_tmp = s_ip_coordinates;
71  auto d_quadCoordinates_tmp = d_quadCoordinates;
72  if constexpr(Sacado::IsADType<typename EvalT::ScalarT>::value) {
73  d_quadCoordinates_tmp(i,j,k).val() = s_ip_coordinates_tmp(i,j,k);
74  }
75  else {
76  d_quadCoordinates_tmp(i,j,k) = s_ip_coordinates_tmp(i,j,k);
77  }
78  });
79 }
80 
81 // **********************************************************************
82 #endif
std::vector< int >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
KOKKOS_FORCEINLINE_FUNCTION array_type get_static_view()
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)