Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ResponseScatterEvaluator_IPCoordinates_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_RESPONSE_SCATTER_EVALUATOR_IPCoordinates_IMPL_HPP
12 #define PANZER_RESPONSE_SCATTER_EVALUATOR_IPCoordinates_IMPL_HPP
13 
14 #include <iostream>
15 #include <string>
16 
17 #include "PanzerDiscFE_config.hpp"
18 
19 #include "Phalanx_Evaluator_Macros.hpp"
20 #include "Phalanx_MDField.hpp"
21 #include "Phalanx_DataLayout_MDALayout.hpp"
22 
23 #include "Panzer_ResponseBase.hpp"
24 #include "Panzer_Dimension.hpp"
27 
28 namespace panzer {
29 
33 template<typename EvalT, typename Traits>
35 ResponseScatterEvaluator_IPCoordinates(const std::string & name,
36  int ir_order)
37  : responseName_(name), ir_order_(ir_order)
38 {
39  using Teuchos::RCP;
40  using Teuchos::rcp;
41 
42  std::string dummyName = ResponseBase::buildLookupName(name) + " dummy target";
43 
44  // build dummy target tag
46  scatterHolder_ = rcp(new PHX::Tag<ScalarT>(dummyName,dl_dummy));
47  this->addEvaluatedField(*scatterHolder_);
48 
49  std::string n = "IPCoordinates Response Scatter: " + name;
50  this->setName(n);
51 }
52 
53 template<typename EvalT, typename Traits>
56 {
57  // extract response object
58  responseObj_ = Teuchos::rcp_dynamic_cast<Response_IPCoordinates<EvalT> >(
59  d.gedc->getDataObject(ResponseBase::buildLookupName(responseName_)),true);
60 }
61 
62 
63 template<typename EvalT, typename Traits>
66  PHX::FieldManager<Traits>& /* fm */)
67 {
68  ir_index_ = panzer::getIntegrationRuleIndex(ir_order_,(*sd.worksets_)[0], this->wda);
69 }
70 
71 template<typename EvalT, typename Traits>
74 {
75  // Kokkos::DynRankView<double,PHX::Device>& workset_coords = (this->wda(workset).int_rules[ir_index_])->ip_coordinates;
76  IntegrationValues2<double> & iv = *this->wda(workset).int_rules[ir_index_];
77 
78  if (tmpCoords_.size() != Teuchos::as<std::size_t>(iv.ip_coordinates.extent(2))) {
79  tmpCoords_.resize(iv.ip_coordinates.extent(2));
80  for(std::size_t dim=0;dim<tmpCoords_.size();dim++)
81  tmpCoords_[dim].clear();
82  }
83 
84  auto ip_coordinates_h = Kokkos::create_mirror_view(PHX::as_view(iv.ip_coordinates));
85  Kokkos::deep_copy(ip_coordinates_h, PHX::as_view(iv.ip_coordinates));
86 
87  // This ordering is for the DataTransferKit. It blocks all x
88  // coordinates for a set of points, then all y coordinates and if
89  // required all z coordinates.
90  for (int dim = 0; dim < iv.ip_coordinates.extent_int(2); ++dim)
91  for (index_t cell = 0; cell < workset.num_cells; ++cell)
92  for (int ip = 0; ip < iv.ip_coordinates.extent_int(1); ++ip)
93  tmpCoords_[dim].push_back(ip_coordinates_h(static_cast<int>(cell),ip,dim));
94 }
95 
96 //**********************************************************************
97 template<typename EvalT, typename Traits>
99 postEvaluate(typename Traits::PostEvalData /* data */)
100 {
101  std::vector<panzer::Traits::Residual::ScalarT> & coords = *responseObj_->getNonconstCoords();
102  coords.clear();
103 
104  for (std::size_t dim = 0; dim < tmpCoords_.size(); ++dim) {
105  for (typename std::vector<ScalarT>::const_iterator x=tmpCoords_[dim].begin(); x != tmpCoords_[dim].end(); ++ x)
106  coords.push_back(Sacado::scalarValue(*x));
107  }
108 
109  tmpCoords_.clear();
110 }
111 
112 }
113 
114 #endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Teuchos::RCP< GlobalEvaluationDataContainer > gedc
int num_cells
DEPRECATED - use: numCells()
static std::string buildLookupName(const std::string &responseName)
ResponseScatterEvaluator_IPCoordinates(const std::string &name, int ir_order)
A constructor with concrete arguments instead of a parameter list.
std::vector< int >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_