11 #ifndef __Panzer_ReorderADValues_Evaluator_impl_hpp__
12 #define __Panzer_ReorderADValues_Evaluator_impl_hpp__
16 #include "Teuchos_Assert.hpp"
17 #include "Teuchos_FancyOStream.hpp"
21 #include "Phalanx_DataLayout.hpp"
23 template<
typename EvalT,
typename TRAITS>
26 const std::vector<std::string> & inFieldNames,
35 for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
40 this->addDependentField(inFields_[eq]);
41 this->addEvaluatedField(outFields_[eq]);
44 this->setName(outPrefix+
" Reorder AD Values");
48 template<
typename EvalT,
typename TRAITS>
51 const std::vector<std::string> & inFieldNames,
52 const std::vector<std::string> & inDOFs,
53 const std::vector<std::string> & outDOFs,
63 for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
68 this->addDependentField(inFields_[eq]);
69 this->addEvaluatedField(outFields_[eq]);
72 this->setName(
"Reorder AD Values");
76 template<
typename EvalT,
typename TRAITS>
80 for(std::size_t i = 0; i < inFields_.size(); ++i)
81 outFields_[i].deep_copy(inFields_[i]);
88 template<
typename TRAITS>
91 const std::vector<std::string> & inFieldNames,
93 const std::string & elementBlock,
102 for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
108 this->addEvaluatedField(outFields_[eq]);
111 buildSrcToDestMap(elementBlock,
115 this->setName(outPrefix+
" Reorder AD Values");
120 template<
typename TRAITS>
123 const std::vector<std::string> & inFieldNames,
124 const std::vector<std::string> & inDOFs,
125 const std::vector<std::string> & outDOFs,
127 const std::string & elementBlock,
135 std::map<int,int> fieldNumberMaps;
136 for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
141 this->addDependentField(inFields_[eq]);
142 this->addEvaluatedField(outFields_[eq]);
144 this->addUnsharedField(outFields_[eq].fieldTag().clone());
148 for(std::size_t i=0;i<inDOFs.size();i++) {
149 int srcFieldNum = indexerSrc.
getFieldNum(inDOFs[i]);
150 int dstFieldNum = indexerDest.
getFieldNum(outDOFs[i]);
154 fieldNumberMaps[srcFieldNum] = dstFieldNum;
157 buildSrcToDestMap(elementBlock,
162 this->setName(
"Reorder AD Values");
166 template<
typename TRAITS>
171 for(std::size_t fieldIndex = 0; fieldIndex < inFields_.size(); ++fieldIndex) {
173 const auto & inField_v = inFields_[fieldIndex].get_view();
174 const auto & outField_v = outFields_[fieldIndex].get_view();
175 const auto & dstFromSrcMap_v = dstFromSrcMapView_;
177 if(inField_v.size()>0) {
179 const auto rank = inField_v.rank();
181 Kokkos::parallel_for(
"ReorderADValues: Jacobian rank 2",workset.num_cells,KOKKOS_LAMBDA(
const int& i){
182 outField_v(i).val() = inField_v(i).val();
183 for (
size_t dx = 0; dx < dstFromSrcMap_v.size(); ++dx)
184 outField_v(i).fastAccessDx(dx) = inField_v(i).fastAccessDx(dstFromSrcMap_v(dx));
188 Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<2>> policy({0,0},{
static_cast<int64_t
>(workset.num_cells),static_cast<int64_t>(inField_v.extent(1))});
189 Kokkos::parallel_for(
"ReorderADValues: Jacobian rank 2",policy,KOKKOS_LAMBDA(
const int& i,
const int& j){
190 outField_v(i,j).val() = inField_v(i,j).val();
191 for (
size_t dx = 0; dx < dstFromSrcMap_v.size(); ++dx)
192 outField_v(i,j).fastAccessDx(dx) = inField_v(i,j).fastAccessDx(dstFromSrcMap_v(dx));
196 Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<3>> policy({0,0,0},{
static_cast<int64_t
>(workset.num_cells),
197 static_cast<int64_t>(inField_v.extent(1)),static_cast<int64_t>(inField_v.extent(2))});
198 Kokkos::parallel_for(
"ReorderADValues: Jacobian rank 2",policy,KOKKOS_LAMBDA(
const int& i,
const int& j,
const int& k){
199 outField_v(i,j,k).val() = inField_v(i,j,k).val();
200 for (
size_t dx = 0; dx < dstFromSrcMap_v.size(); ++dx)
201 outField_v(i,j,k).fastAccessDx(dx) = inField_v(i,j,k).fastAccessDx(dstFromSrcMap_v(dx));
205 Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<4>> policy({0,0,0,0},{
static_cast<int64_t
>(workset.num_cells),
206 static_cast<int64_t>(inField_v.extent(1)),static_cast<int64_t>(inField_v.extent(2)),
207 static_cast<int64_t>(inField_v.extent(3))});
208 Kokkos::parallel_for(
"ReorderADValues: Jacobian rank 2",policy,KOKKOS_LAMBDA(
const int& i,
const int& j,
const int& k,
const int& l){
209 outField_v(i,j,k,l).val() = inField_v(i,j,k,l).val();
210 for (
size_t dx = 0; dx < dstFromSrcMap_v.size(); ++dx)
211 outField_v(i,j,k,l).fastAccessDx(dx) = inField_v(i,j,k,l).fastAccessDx(dstFromSrcMap_v(dx));
215 Kokkos::MDRangePolicy<PHX::Device,Kokkos::Rank<5>> policy({0,0,0,0,0},{
static_cast<int64_t
>(workset.num_cells),
216 static_cast<int64_t>(inField_v.extent(1)),static_cast<int64_t>(inField_v.extent(2)),
217 static_cast<int64_t>(inField_v.extent(3)),static_cast<int64_t>(inField_v.extent(4))});
218 Kokkos::parallel_for(
"ReorderADValues: Jacobian rank 2",policy,KOKKOS_LAMBDA(
const int& i,
const int& j,
const int& k,
const int& l,
const int& m){
219 outField_v(i,j,k,l,m).val() = inField_v(i,j,k,l,m).val();
220 for (
size_t dx = 0; dx < dstFromSrcMap_v.size(); ++dx)
221 outField_v(i,j,k,l,m).fastAccessDx(dx) = inField_v(i,j,k,l,m).fastAccessDx(dstFromSrcMap_v(dx));
232 template<
typename TRAITS>
247 std::map<int,int> fieldNumberMaps;
248 for(std::size_t i=0;i<dstFieldsNum.size();i++) {
249 std::string fieldName = indexerDest.
getFieldString(dstFieldsNum[i]);
251 int srcFieldNum = indexerSrc.
getFieldNum(fieldName);
253 fieldNumberMaps[srcFieldNum] = dstFieldsNum[i];
255 out <<
"Warning: Reorder AD Values can't find field \"" << fieldName <<
"\"" << std::endl;
258 buildSrcToDestMap(elementBlock,fieldNumberMaps,indexerSrc,indexerDest);
262 template<
typename TRAITS>
265 const std::map<int,int> & fieldNumberMaps,
270 std::map<int,int> offsetMap;
271 for(std::map<int,int>::const_iterator itr=fieldNumberMaps.begin();
272 itr!=fieldNumberMaps.end();++itr) {
273 int srcField = itr->first;
274 int dstField = itr->second;
276 const std::vector<int> & srcOffsets = indexerSrc.
getGIDFieldOffsets(elementBlock,srcField);
277 const std::vector<int> & dstOffsets = indexerDest.
getGIDFieldOffsets(elementBlock,dstField);
281 for(std::size_t i=0;i<srcOffsets.size();i++) {
282 offsetMap[srcOffsets[i]] = dstOffsets[i];
286 maxDest = dstOffsets[i]>maxDest ? dstOffsets[i] : maxDest;
292 std::vector<int> dstFromSrcMap(maxDest+1,-1);
293 for(std::map<int,int>::const_iterator itr=offsetMap.begin();
294 itr!=offsetMap.end();++itr) {
295 dstFromSrcMap[itr->second] = itr->first;
299 auto dfsm_host = Kokkos::create_mirror_view(Kokkos::HostSpace(),dstFromSrcMapView_);
300 for (
size_t i=0; i < dstFromSrcMapView_.size(); ++i)
301 dfsm_host(i) = dstFromSrcMap[i];
303 Kokkos::deep_copy(dstFromSrcMapView_,dfsm_host);
virtual const std::string & getFieldString(int num) const =0
Reverse lookup of the field string from a field number.
virtual Teuchos::RCP< Teuchos::Comm< int > > getComm() const =0
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::vector< PHX::MDField< ScalarT > > outFields_
void evaluateFields(typename TRAITS::EvalData d)
std::vector< PHX::MDField< const ScalarT > > inFields_
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
Reorders the ad values of a specified field to match a different unique global indexer.
virtual const std::vector< int > & getGIDFieldOffsets(const std::string &blockId, int fieldNum) const =0
Use the field pattern so that you can find a particular field in the GIDs array.
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0
ReorderADValues_Evaluator(const std::string &outPrefix, const std::vector< std::string > &inFieldNames, const std::vector< Teuchos::RCP< PHX::DataLayout > > &fieldLayouts, const std::string &elementBlock, const GlobalIndexer &indexerSrc, const GlobalIndexer &indexerDest)
#define TEUCHOS_ASSERT(assertion_test)