Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ReorderADValues_Evaluator_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_ReorderADValues_Evaluator_impl_hpp__
12 #define __Panzer_ReorderADValues_Evaluator_impl_hpp__
13 
14 
15 #include "Teuchos_RCP.hpp"
16 #include "Teuchos_Assert.hpp"
17 #include "Teuchos_FancyOStream.hpp"
18 
19 #include "Panzer_GlobalIndexer.hpp"
20 
21 #include "Phalanx_DataLayout.hpp"
22 
23 template<typename EvalT,typename TRAITS>
25 ReorderADValues_Evaluator(const std::string & outPrefix,
26  const std::vector<std::string> & inFieldNames,
27  const std::vector<Teuchos::RCP<PHX::DataLayout> > & fieldLayouts,
28  const std::string & /* elementBlock */,
29  const GlobalIndexer & /* indexerSrc */,
30  const GlobalIndexer & /* indexerDest */)
31 {
32  TEUCHOS_ASSERT(inFieldNames.size()==fieldLayouts.size());
33 
34  // build the vector of fields that this is dependent on
35  for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
36  inFields_.push_back(PHX::MDField<const ScalarT>(inFieldNames[eq],fieldLayouts[eq]));
37  outFields_.push_back(PHX::MDField<ScalarT>(outPrefix+inFieldNames[eq],fieldLayouts[eq]));
38 
39  // tell the field manager that we depend on this field
40  this->addDependentField(inFields_[eq]);
41  this->addEvaluatedField(outFields_[eq]);
42  }
43 
44  this->setName(outPrefix+" Reorder AD Values");
45 }
46 
47 // **********************************************************************
48 template<typename EvalT,typename TRAITS>
50 ReorderADValues_Evaluator(const std::string & outPrefix,
51  const std::vector<std::string> & inFieldNames,
52  const std::vector<std::string> & inDOFs,
53  const std::vector<std::string> & outDOFs,
54  const std::vector<Teuchos::RCP<PHX::DataLayout> > & fieldLayouts,
55  const std::string & /* elementBlock */,
56  const GlobalIndexer & /* indexerSrc */,
57  const GlobalIndexer & /* indexerDest */)
58 {
59  TEUCHOS_ASSERT(inFieldNames.size()==fieldLayouts.size());
60  TEUCHOS_ASSERT(inDOFs.size()==outDOFs.size());
61 
62  // build the vector of fields that this is dependent on
63  for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
64  inFields_.push_back(PHX::MDField<const ScalarT>(inFieldNames[eq],fieldLayouts[eq]));
65  outFields_.push_back(PHX::MDField<ScalarT>(outPrefix+inFieldNames[eq],fieldLayouts[eq]));
66 
67  // tell the field manager that we depend on this field
68  this->addDependentField(inFields_[eq]);
69  this->addEvaluatedField(outFields_[eq]);
70  }
71 
72  this->setName("Reorder AD Values");
73 }
74 
75 // **********************************************************************
76 template<typename EvalT,typename TRAITS>
78 evaluateFields(typename TRAITS::EvalData /* workset */)
79 {
80  for(std::size_t i = 0; i < inFields_.size(); ++i)
81  outFields_[i].deep_copy(inFields_[i]);
82 }
83 
84 // **********************************************************************
85 // Jacobian
86 // **********************************************************************
87 
88 template<typename TRAITS>
90 ReorderADValues_Evaluator(const std::string & outPrefix,
91  const std::vector<std::string> & inFieldNames,
92  const std::vector<Teuchos::RCP<PHX::DataLayout> > & fieldLayouts,
93  const std::string & elementBlock,
94  const GlobalIndexer & indexerSrc,
95  const GlobalIndexer & indexerDest)
96 {
97  TEUCHOS_ASSERT(inFieldNames.size()==fieldLayouts.size());
98 
99  // build the vector of fields that this is dependent on
100  inFields_.resize(inFieldNames.size());
101  outFields_.resize(inFieldNames.size());
102  for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
103  inFields_[eq] = PHX::MDField<const ScalarT>(inFieldNames[eq],fieldLayouts[eq]);
104  outFields_[eq] = PHX::MDField<ScalarT>(outPrefix+inFieldNames[eq],fieldLayouts[eq]);
105 
106  // tell the field manager that we depend on this field
107  this->addDependentField(inFields_[eq]);
108  this->addEvaluatedField(outFields_[eq]);
109  }
110 
111  buildSrcToDestMap(elementBlock,
112  indexerSrc,
113  indexerDest);
114 
115  this->setName(outPrefix+" Reorder AD Values");
116 }
117 
118 // **********************************************************************
119 
120 template<typename TRAITS>
122 ReorderADValues_Evaluator(const std::string & outPrefix,
123  const std::vector<std::string> & inFieldNames,
124  const std::vector<std::string> & inDOFs,
125  const std::vector<std::string> & outDOFs,
126  const std::vector<Teuchos::RCP<PHX::DataLayout> > & fieldLayouts,
127  const std::string & elementBlock,
128  const GlobalIndexer & indexerSrc,
129  const GlobalIndexer & indexerDest)
130 {
131  TEUCHOS_ASSERT(inFieldNames.size()==fieldLayouts.size());
132  TEUCHOS_ASSERT(inDOFs.size()==outDOFs.size());
133 
134  // build the vector of fields that this is dependent on
135  std::map<int,int> fieldNumberMaps;
136  for (std::size_t eq = 0; eq < inFieldNames.size(); ++eq) {
137  inFields_.push_back(PHX::MDField<const ScalarT>(inFieldNames[eq],fieldLayouts[eq]));
138  outFields_.push_back(PHX::MDField<ScalarT>(outPrefix+inFieldNames[eq],fieldLayouts[eq]));
139 
140  // tell the field manager that we depend on this field
141  this->addDependentField(inFields_[eq]);
142  this->addEvaluatedField(outFields_[eq]);
143  // Don't share so we can avoid zeroing out off blck Jacobian entries
144  this->addUnsharedField(outFields_[eq].fieldTag().clone());
145  }
146 
147  // build a int-int map that associates fields
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]);
151  TEUCHOS_ASSERT(srcFieldNum>=0);
152  TEUCHOS_ASSERT(dstFieldNum>=0);
153 
154  fieldNumberMaps[srcFieldNum] = dstFieldNum;
155  }
156 
157  buildSrcToDestMap(elementBlock,
158  fieldNumberMaps,
159  indexerSrc,
160  indexerDest);
161 
162  this->setName("Reorder AD Values");
163 }
164 
165 // **********************************************************************
166 template<typename TRAITS>
168 evaluateFields(typename TRAITS::EvalData workset)
169 {
170  // for AD data do a reordering
171  for(std::size_t fieldIndex = 0; fieldIndex < inFields_.size(); ++fieldIndex) {
172 
173  const auto & inField_v = inFields_[fieldIndex].get_view();
174  const auto & outField_v = outFields_[fieldIndex].get_view();
175  const auto & dstFromSrcMap_v = dstFromSrcMapView_;
176 
177  if(inField_v.size()>0) {
178 
179  const auto rank = inField_v.rank();
180  if (rank==1) {
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));
185  });
186  }
187  else if (rank==2) {
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));
193  });
194  }
195  else if (rank==3) {
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));
202  });
203  }
204  else if (rank==4) {
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));
212  });
213  }
214  else if (rank==5) {
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));
222  });
223  }
224  else {
225  TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"ERROR AD Reorder, rank size " << rank << " not supported!");
226  }
227  }
228  }
229 }
230 
231 // **********************************************************************
232 template<typename TRAITS>
234 buildSrcToDestMap(const std::string & elementBlock,
235  const GlobalIndexer & indexerSrc,
236  const GlobalIndexer & indexerDest)
237 {
238  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
239  out.setOutputToRootOnly(0);
240 
241  TEUCHOS_ASSERT(indexerSrc.getComm()!=Teuchos::null);
242  TEUCHOS_ASSERT(indexerDest.getComm()!=Teuchos::null);
243 
244  const std::vector<int> & dstFieldsNum = indexerDest.getBlockFieldNumbers(elementBlock);
245 
246  // build a map between destination field numbers and source field numbers
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]);
250 
251  int srcFieldNum = indexerSrc.getFieldNum(fieldName);
252  if(srcFieldNum>=0)
253  fieldNumberMaps[srcFieldNum] = dstFieldsNum[i];
254  else
255  out << "Warning: Reorder AD Values can't find field \"" << fieldName << "\"" << std::endl;
256  }
257 
258  buildSrcToDestMap(elementBlock,fieldNumberMaps,indexerSrc,indexerDest);
259 }
260 
261 // **********************************************************************
262 template<typename TRAITS>
264 buildSrcToDestMap(const std::string & elementBlock,
265  const std::map<int,int> & fieldNumberMaps,
266  const GlobalIndexer & indexerSrc,
267  const GlobalIndexer & indexerDest)
268 {
269  int maxDest = -1;
270  std::map<int,int> offsetMap; // map from source to destination offsets
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;
275 
276  const std::vector<int> & srcOffsets = indexerSrc.getGIDFieldOffsets(elementBlock,srcField);
277  const std::vector<int> & dstOffsets = indexerDest.getGIDFieldOffsets(elementBlock,dstField);
278 
279  // field should be the same size
280  TEUCHOS_ASSERT(srcOffsets.size()==dstOffsets.size());
281  for(std::size_t i=0;i<srcOffsets.size();i++) {
282  offsetMap[srcOffsets[i]] = dstOffsets[i];
283 
284  // provides a size for allocating an array below: we will be able
285  // to index into dstFromSrcMap_ in a simple way
286  maxDest = dstOffsets[i]>maxDest ? dstOffsets[i] : maxDest;
287  }
288  }
289 
290  // Build map
291  TEUCHOS_ASSERT(maxDest>0);
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;
296  }
297 
298  dstFromSrcMapView_ = Kokkos::View<int*>("dstFromSrcMapView_",dstFromSrcMap.size());
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];
302 
303  Kokkos::deep_copy(dstFromSrcMapView_,dfsm_host);
304 }
305 
306 // **********************************************************************
307 
308 #endif
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_
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)