Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_L2Projection.cpp
Go to the documentation of this file.
1 // @HEADER
2 // @HEADER
3 
4 #ifndef PANZER_L2_PROJECTION_IMPL_HPP
5 #define PANZER_L2_PROJECTION_IMPL_HPP
6 
8 #include "Tpetra_CrsGraph.hpp"
9 #include "Tpetra_MultiVector.hpp"
10 #include "Tpetra_CrsMatrix.hpp"
11 #include "Panzer_L2Projection.hpp"
14 #include "Panzer_TpetraLinearObjFactory.hpp"
22 #include "Panzer_Workset.hpp"
23 
24 namespace panzer {
25 
27  const panzer::IntegrationDescriptor& integrationDescriptor,
28  const Teuchos::RCP<const Teuchos::MpiComm<int>>& comm,
29  const Teuchos::RCP<const panzer::ConnManager>& connManager,
30  const std::vector<std::string>& elementBlockNames,
31  const Teuchos::RCP<panzer::WorksetContainer> worksetContainer)
32  {
33  targetBasisDescriptor_ = targetBasis;
34  integrationDescriptor_ = integrationDescriptor;
35  comm_ = comm;
36  connManager_ = connManager;
37  elementBlockNames_ = elementBlockNames;
38  worksetContainer_ = worksetContainer;
39  setupCalled_ = true;
40 
41  // Build target DOF Manager
42  targetGlobalIndexer_ =
43  Teuchos::rcp(new panzer::DOFManager(Teuchos::rcp_const_cast<panzer::ConnManager>(connManager),*(comm->getRawMpiComm())));
44 
45  // For hybrid mesh, blocks could have different topology
46  for (const auto& eBlock : elementBlockNames_) {
47  std::vector<shards::CellTopology> topologies;
48  connManager_->getElementBlockTopologies(topologies);
49  std::vector<std::string> ebNames;
50  connManager_->getElementBlockIds(ebNames);
51  const auto search = std::find(ebNames.cbegin(),ebNames.cend(),eBlock);
52  TEUCHOS_ASSERT(search != ebNames.cend());
53  const int index = std::distance(ebNames.cbegin(),search);
54  const auto& cellTopology = topologies[index];
55 
56  auto intrepidBasis = panzer::createIntrepid2Basis<PHX::Device,double,double>(targetBasisDescriptor_.getType(),
57  targetBasisDescriptor_.getOrder(),
58  cellTopology);
60  // Field name is the basis type
61  targetGlobalIndexer_->addField(eBlock,targetBasisDescriptor_.getType(),fieldPattern);
62  }
63 
64  targetGlobalIndexer_->buildGlobalUnknowns();
65 
66  // Check workset needs are correct
67  }
68 
71  {return targetGlobalIndexer_;}
72 
75  const std::unordered_map<std::string,double>* elementBlockMultipliers)
76  {
77  PANZER_FUNC_TIME_MONITOR_DIFF("L2Projection::Build Mass Matrix",TopBuildMassMatrix);
78  TEUCHOS_ASSERT(setupCalled_);
79 
80  if (elementBlockMultipliers != nullptr) {
81  TEUCHOS_ASSERT(elementBlockMultipliers->size() == elementBlockNames_.size());
82  }
83 
84  // Allocate the owned matrix
85  std::vector<Teuchos::RCP<const panzer::GlobalIndexer>> indexers;
86  indexers.push_back(targetGlobalIndexer_);
87 
89 
90  auto ownedMatrix = factory.getTpetraMatrix(0,0);
91  auto ghostedMatrix = factory.getGhostedTpetraMatrix(0,0);
92  ownedMatrix->resumeFill();
93  ownedMatrix->setAllToScalar(0.0);
94  ghostedMatrix->resumeFill();
95  ghostedMatrix->setAllToScalar(0.0);
96  PHX::Device().fence();
97 
98  auto M = ghostedMatrix->getLocalMatrix();
99  const int fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
100 
101  const bool is_scalar = targetBasisDescriptor_.getType()=="HGrad" || targetBasisDescriptor_.getType()=="Const" || targetBasisDescriptor_.getType()=="HVol";
102 
103  // Loop over element blocks and fill mass matrix
104  if(is_scalar){
105  for (const auto& block : elementBlockNames_) {
106 
107  double ebMultiplier = 1.0;
108  if (elementBlockMultipliers != nullptr)
109  ebMultiplier = elementBlockMultipliers->find(block)->second;
110 
111  // Based on descriptor, currently assumes there should only be one workset
113  const auto worksets = worksetContainer_->getWorksets(wd);
114 
115  for (const auto& workset : *worksets) {
116 
117  const auto basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
118 
119  const auto unweightedBasis = basisValues.basis_scalar;
120  const auto weightedBasis = basisValues.weighted_basis_scalar;
121 
122  // Offsets (this assumes UVM, need to fix)
123  const std::vector<panzer::LocalOrdinal>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
124  PHX::View<panzer::LocalOrdinal*> kOffsets("MassMatrix: Offsets",offsets.size());
125  for(const auto& i : offsets)
126  kOffsets(i) = offsets[i];
127 
128  PHX::Device().fence();
129 
130  // Local Ids
131  Kokkos::View<panzer::LocalOrdinal**,PHX::Device> localIds("MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
132  targetGlobalIndexer_->getElementBlockGIDCount(block));
133 
134  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
135  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
136 
137  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
138 
139  const int numBasisPoints = static_cast<int>(weightedBasis.extent(1));
140  if ( use_lumping ) {
141  Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (const int& cell) {
142  double total_mass = 0.0, trace = 0.0;
143 
144  panzer::LocalOrdinal cLIDs[256];
145  const int numIds = static_cast<int>(localIds.extent(1));
146  for(int i=0;i<numIds;++i)
147  cLIDs[i] = localIds(cell,i);
148 
149  double vals[256]={0.0};
150  const int numQP = static_cast<int>(unweightedBasis.extent(2));
151 
152  for (int row=0; row < numBasisPoints; ++row) {
153  for (int col=0; col < numIds; ++col) {
154  for (int qp=0; qp < numQP; ++qp) {
155  auto tmp = unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier;
156  total_mass += tmp;
157  if (col == row )
158  trace += tmp;
159  }
160  }
161  }
162 
163  for (int row=0; row < numBasisPoints; ++row) {
164  for (int col=0; col < numBasisPoints; ++col)
165  vals[col] = 0;
166 
167  int offset = kOffsets(row);
168  panzer::LocalOrdinal lid = localIds(cell,offset);
169  int col = row;
170  vals[col] = 0.0;
171  for (int qp=0; qp < numQP; ++qp)
172  vals[col] += unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier * total_mass / trace;
173 
174  M.sumIntoValues(lid,cLIDs,numIds,vals,true,true);
175  }
176  });
177 
178  } else {
179  Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (const int& cell) {
180 
181  panzer::LocalOrdinal cLIDs[256];
182  const int numIds = static_cast<int>(localIds.extent(1));
183  for(int i=0;i<numIds;++i)
184  cLIDs[i] = localIds(cell,i);
185 
186  double vals[256];
187  const int numQP = static_cast<int>(unweightedBasis.extent(2));
188 
189  for (int row=0; row < numBasisPoints; ++row) {
190  int offset = kOffsets(row);
191  panzer::LocalOrdinal lid = localIds(cell,offset);
192 
193  for (int col=0; col < numIds; ++col) {
194  vals[col] = 0.0;
195  for (int qp=0; qp < numQP; ++qp)
196  vals[col] += unweightedBasis(cell,row,qp) * weightedBasis(cell,col,qp) * ebMultiplier;
197  }
198  M.sumIntoValues(lid,cLIDs,numIds,vals,true,true);
199 
200  }
201 
202  });
203  }
204  }
205  }
206  } else {
207  for (const auto& block : elementBlockNames_) {
208 
209  double ebMultiplier = 1.0;
210  if (elementBlockMultipliers != nullptr)
211  ebMultiplier = elementBlockMultipliers->find(block)->second;
212 
213  // Based on descriptor, currently assumes there should only be one workset
215  const auto& worksets = worksetContainer_->getWorksets(wd);
216 
217  for (const auto& workset : *worksets) {
218 
219  const auto basisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
220 
221  const auto unweightedBasis = basisValues.basis_vector;
222  const auto weightedBasis = basisValues.weighted_basis_vector;
223 
224  // Offsets (this assumes UVM, need to fix)
225  const std::vector<panzer::LocalOrdinal>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
226  PHX::View<panzer::LocalOrdinal*> kOffsets("MassMatrix: Offsets",offsets.size());
227  for(const auto& i : offsets)
228  kOffsets(i) = offsets[i];
229 
230  PHX::Device().fence();
231 
232  // Local Ids
233  Kokkos::View<panzer::LocalOrdinal**,PHX::Device> localIds("MassMatrix: LocalIds", workset.numOwnedCells()+workset.numGhostCells()+workset.numVirtualCells(),
234  targetGlobalIndexer_->getElementBlockGIDCount(block));
235 
236  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
237  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
238 
239  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,localIds);
240 
241  const int numBasisPoints = static_cast<int>(weightedBasis.extent(1));
242  Kokkos::parallel_for(workset.numOwnedCells(),KOKKOS_LAMBDA (const int& cell) {
243 
244  panzer::LocalOrdinal cLIDs[256];
245  const int numIds = static_cast<int>(localIds.extent(1));
246  for(int i=0;i<numIds;++i)
247  cLIDs[i] = localIds(cell,i);
248 
249  double vals[256];
250  const int numQP = static_cast<int>(unweightedBasis.extent(2));
251 
252  for (int qp=0; qp < numQP; ++qp) {
253  for (int row=0; row < numBasisPoints; ++row) {
254  int offset = kOffsets(row);
255  panzer::LocalOrdinal lid = localIds(cell,offset);
256 
257  for (int col=0; col < numIds; ++col){
258  vals[col] = 0.0;
259  for(int dim=0; dim < static_cast<int>(weightedBasis.extent(3)); ++dim)
260  vals[col] += unweightedBasis(cell,row,qp,dim) * weightedBasis(cell,col,qp,dim) * ebMultiplier;
261  }
262 
263  M.sumIntoValues(lid,cLIDs,numIds,vals,true,true);
264  }
265  }
266 
267  });
268  }
269  }
270  }
271  PHX::exec_space().fence();
272 
273  {
274  PANZER_FUNC_TIME_MONITOR_DIFF("Exporting of mass matrix",ExportMM);
275  auto map = factory.getMap(0);
276  ghostedMatrix->fillComplete(map,map);
277  const auto exporter = factory.getGhostedExport(0);
278  ownedMatrix->doExport(*ghostedMatrix, *exporter, Tpetra::ADD);
279  ownedMatrix->fillComplete();
280  }
281  return ownedMatrix;
282  }
283 
286  {
287  PANZER_FUNC_TIME_MONITOR_DIFF("L2Projection<panzer::LocalOrdinal,panzer::GlobalOrdinal>::buildInverseLumpedMassMatrix",buildInvLMM);
288  using Teuchos::rcp;
289  const auto massMatrix = this->buildMassMatrix(true);
290  const auto lumpedMassMatrix = rcp(new Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>(massMatrix->getDomainMap(),1,true));
291  const auto tmp = rcp(new Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>(massMatrix->getRangeMap(),1,false));
292  tmp->putScalar(1.0);
293  {
294  PANZER_FUNC_TIME_MONITOR_DIFF("Apply",Apply);
295  massMatrix->apply(*tmp,*lumpedMassMatrix);
296  }
297  {
298  PANZER_FUNC_TIME_MONITOR_DIFF("Reciprocal",Reciprocal);
299  lumpedMassMatrix->reciprocal(*lumpedMassMatrix);
300  }
301  return lumpedMassMatrix;
302  }
303 
306  const Teuchos::RCP<const Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>>& inputOwnedSourceMap,
307  const std::string& sourceFieldName,
308  const panzer::BasisDescriptor& sourceBasisDescriptor,
309  const int directionIndex)
310  {
311  // *******************
312  // Create Retangular matrix (both ghosted and owned).
313  // *******************
314  using Teuchos::RCP;
315  using Teuchos::rcp;
316  using MapType = Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
317  using GraphType = Tpetra::CrsGraph<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
318  using ExportType = Tpetra::Export<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
319  using MatrixType = Tpetra::CrsMatrix<double,panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
320 
321  // *******************
322  // Build ghosted graph
323  // *******************
324 
325  RCP<MapType> ghostedTargetMap;
326  {
327  std::vector<panzer::GlobalOrdinal> indices;
328  targetGlobalIndexer_->getOwnedAndGhostedIndices(indices);
329  ghostedTargetMap = rcp(new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
330  }
331 
332  RCP<MapType> ghostedSourceMap;
333  {
334  std::vector<panzer::GlobalOrdinal> indices;
335  sourceGlobalIndexer.getOwnedAndGhostedIndices(indices);
336  ghostedSourceMap = rcp(new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
337  }
338 
339 
340  // Now insert the non-zero pattern per row
341  // count number of entries per row; required by CrsGraph constructor
342  std::vector<size_t> nEntriesPerRow(ghostedTargetMap->getNodeNumElements(),0);
343  std::vector<std::string> elementBlockIds;
344  targetGlobalIndexer_->getElementBlockIds(elementBlockIds);
345  std::vector<std::string>::const_iterator blockItr;
346  for (blockItr=elementBlockIds.begin();blockItr!=elementBlockIds.end();++blockItr) {
347  std::string blockId = *blockItr;
348  const std::vector<panzer::LocalOrdinal> & elements = targetGlobalIndexer_->getElementBlock(blockId);
349 
350  std::vector<panzer::GlobalOrdinal> row_gids;
351  std::vector<panzer::GlobalOrdinal> col_gids;
352 
353  for(std::size_t elmt=0;elmt<elements.size();elmt++) {
354  targetGlobalIndexer_->getElementGIDs(elements[elmt],row_gids);
355  sourceGlobalIndexer.getElementGIDs(elements[elmt],col_gids);
356  for(std::size_t row=0;row<row_gids.size();row++) {
357  panzer::LocalOrdinal lid =
358  ghostedTargetMap->getLocalElement(row_gids[row]);
359  nEntriesPerRow[lid] += col_gids.size();
360  }
361  }
362  }
363 
364  Teuchos::ArrayView<const size_t> nEntriesPerRowView(nEntriesPerRow);
365  RCP<GraphType> ghostedGraph = rcp(new GraphType(ghostedTargetMap,ghostedSourceMap,nEntriesPerRowView,Tpetra::StaticProfile));
366 
367  for (blockItr=elementBlockIds.begin();blockItr!=elementBlockIds.end();++blockItr) {
368  std::string blockId = *blockItr;
369  const std::vector<panzer::LocalOrdinal> & elements = targetGlobalIndexer_->getElementBlock(blockId);
370 
371  std::vector<panzer::GlobalOrdinal> row_gids;
372  std::vector<panzer::GlobalOrdinal> col_gids;
373 
374  for(std::size_t elmt=0;elmt<elements.size();elmt++) {
375  targetGlobalIndexer_->getElementGIDs(elements[elmt],row_gids);
376  sourceGlobalIndexer.getElementGIDs(elements[elmt],col_gids);
377  for(std::size_t row=0;row<row_gids.size();row++)
378  ghostedGraph->insertGlobalIndices(row_gids[row],col_gids);
379  }
380  }
381 
382  RCP<MapType> ownedTargetMap;
383  {
384  std::vector<panzer::GlobalOrdinal> indices;
385  targetGlobalIndexer_->getOwnedIndices(indices);
386  ownedTargetMap = rcp(new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
387  }
388 
389  RCP<const MapType> ownedSourceMap = inputOwnedSourceMap;
390  if (is_null(ownedSourceMap)) {
391  std::vector<panzer::GlobalOrdinal> indices;
392  sourceGlobalIndexer.getOwnedIndices(indices);
393  ownedSourceMap = rcp(new MapType(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),indices,0,comm_));
394  }
395 
396  // Fill complete with owned range and domain map
397  ghostedGraph->fillComplete(ownedSourceMap,ownedTargetMap);
398 
399  // *****************
400  // Build owned graph
401  // *****************
402 
403  RCP<GraphType> ownedGraph = rcp(new GraphType(ownedTargetMap,0));
404  RCP<const ExportType> exporter = rcp(new ExportType(ghostedTargetMap,ownedTargetMap));
405  ownedGraph->doExport(*ghostedGraph, *exporter, Tpetra::INSERT);
406  ownedGraph->fillComplete(ownedSourceMap,ownedTargetMap);
407 
408  RCP<MatrixType> ghostedMatrix = rcp(new MatrixType(ghostedGraph));
409  RCP<MatrixType> ownedMatrix = rcp(new MatrixType(ownedGraph));
410  // ghostedMatrix.fillComplete();
411  // ghostedMatrix.resumeFill();
412 
413  ghostedMatrix->setAllToScalar(0.0);
414  ownedMatrix->setAllToScalar(0.0);
415  PHX::Device().fence();
416 
417  // *******************
418  // Fill ghosted matrix
419  // *******************
420  for (const auto& block : elementBlockNames_) {
421 
423  const auto& worksets = worksetContainer_->getWorksets(wd);
424  for (const auto& workset : *worksets) {
425 
426  // Get target basis values: current implementation assumes target basis is HGrad
427  const auto& targetBasisValues = workset.getBasisValues(targetBasisDescriptor_,integrationDescriptor_);
428  const auto& targetWeightedBasis = targetBasisValues.weighted_basis_scalar.get_static_view();
429 
430  // Sources can be any basis
431  const auto& sourceBasisValues = workset.getBasisValues(sourceBasisDescriptor,integrationDescriptor_);
432  Kokkos::View<const double***,PHX::Device> sourceUnweightedScalarBasis;
433  Kokkos::View<const double****,PHX::Device> sourceUnweightedVectorBasis;
434  bool useRankThreeBasis = false; // default to gradient or vector basis
435  if ( (sourceBasisDescriptor.getType() == "HGrad") || (sourceBasisDescriptor.getType() == "Const") || (sourceBasisDescriptor.getType() == "HVol") ) {
436  if (directionIndex == -1) { // Project dof value
437  sourceUnweightedScalarBasis = sourceBasisValues.basis_scalar.get_static_view();
438  useRankThreeBasis = true;
439  }
440  else { // Project dof gradient of scalar basis
441  sourceUnweightedVectorBasis = sourceBasisValues.grad_basis.get_static_view();
442  }
443  }
444  else { // Project vector value
445  sourceUnweightedVectorBasis = sourceBasisValues.basis_vector.get_static_view();
446  }
447 
448  // Get the element local ids
449  Kokkos::View<panzer::LocalOrdinal**,PHX::Device> targetLocalIds("buildRHSMatrix: targetLocalIds", workset.numOwnedCells(),
450  targetGlobalIndexer_->getElementBlockGIDCount(block));
451  Kokkos::View<panzer::LocalOrdinal**,PHX::Device> sourceLocalIds("buildRHSMatrix: sourceLocalIds", workset.numOwnedCells(),
452  sourceGlobalIndexer.getElementBlockGIDCount(block));
453  {
454  // Remove the ghosted cell ids or the call to getElementLocalIds will spill array bounds
455  const auto cellLocalIdsNoGhost = Kokkos::subview(workset.cell_local_ids_k,std::make_pair(0,workset.numOwnedCells()));
456  targetGlobalIndexer_->getElementLIDs(cellLocalIdsNoGhost,targetLocalIds);
457  sourceGlobalIndexer.getElementLIDs(cellLocalIdsNoGhost,sourceLocalIds);
458  }
459 
460  // Get the offsets
461  Kokkos::View<panzer::LocalOrdinal*,PHX::Device> targetFieldOffsets;
462  {
463  const auto fieldIndex = targetGlobalIndexer_->getFieldNum(targetBasisDescriptor_.getType());
464  const std::vector<panzer::LocalOrdinal>& offsets = targetGlobalIndexer_->getGIDFieldOffsets(block,fieldIndex);
465  targetFieldOffsets = Kokkos::View<panzer::LocalOrdinal*,PHX::Device>("L2Projection:buildRHS:targetFieldOffsets",offsets.size());
466  const auto hostOffsets = Kokkos::create_mirror_view(targetFieldOffsets);
467  for(size_t i=0; i < offsets.size(); ++i)
468  hostOffsets(i) = offsets[i];
469  Kokkos::deep_copy(targetFieldOffsets,hostOffsets);
470  PHX::Device().fence();
471  }
472 
473  Kokkos::View<panzer::LocalOrdinal*,PHX::Device> sourceFieldOffsets;
474  {
475  const auto fieldIndex = sourceGlobalIndexer.getFieldNum(sourceFieldName);
476  const std::vector<panzer::LocalOrdinal>& offsets = sourceGlobalIndexer.getGIDFieldOffsets(block,fieldIndex);
477  TEUCHOS_TEST_FOR_EXCEPTION(offsets.size() == 0, std::runtime_error,
478  "ERROR: panzer::L2Projection::buildRHSMatrix() - The source field, \""
479  << sourceFieldName << "\", does not exist in element block \""
480  << block << "\"!");
481  sourceFieldOffsets = Kokkos::View<panzer::LocalOrdinal*,PHX::Device>("L2Projection:buildRHS:sourceFieldOffsets",offsets.size());
482  const auto hostOffsets = Kokkos::create_mirror_view(sourceFieldOffsets);
483  for(size_t i=0; i <offsets.size(); ++i)
484  hostOffsets(i) = offsets[i];
485  Kokkos::deep_copy(sourceFieldOffsets,hostOffsets);
486  PHX::Device().fence();
487  }
488 
489  const auto localMatrix = ghostedMatrix->getLocalMatrix();
490  const int numRows = static_cast<int>(targetWeightedBasis.extent(1));
491  int tmpNumCols = -1;
492  int tmpNumQP = -1;
493  if (useRankThreeBasis) {
494  tmpNumCols = static_cast<int>(sourceUnweightedScalarBasis.extent(1));
495  tmpNumQP = static_cast<int>(sourceUnweightedScalarBasis.extent(2));
496  }
497  else {
498  tmpNumCols = static_cast<int>(sourceUnweightedVectorBasis.extent(1));
499  tmpNumQP = static_cast<int>(sourceUnweightedVectorBasis.extent(2));
500  }
501  const int numCols = tmpNumCols;
502  const int numQP = tmpNumQP;
503 
504  Kokkos::parallel_for(Kokkos::RangePolicy<PHX::Device>(0,workset.numOwnedCells()),KOKKOS_LAMBDA (const int& cell) {
505  panzer::LocalOrdinal cLIDs[256];
506  double vals[256];
507  for (int row = 0; row < numRows; ++row) {
508  const int rowOffset = targetFieldOffsets(row);
509  const int rowLID = targetLocalIds(cell,rowOffset);
510  for (int col = 0; col < numCols; ++col)
511  vals[col] = 0.0;
512 
513  for (int col = 0; col < numCols; ++col) {
514  for (int qp = 0; qp < numQP; ++qp) {
515  const int colOffset = sourceFieldOffsets(col);
516  const int colLID = sourceLocalIds(cell,colOffset);
517  cLIDs[col] = colLID;
518  if (useRankThreeBasis)
519  vals[col] += sourceUnweightedScalarBasis(cell,col,qp) * targetWeightedBasis(cell,row,qp);
520  else
521  vals[col] += sourceUnweightedVectorBasis(cell,col,qp,directionIndex) * targetWeightedBasis(cell,row,qp);
522  }
523  }
524  localMatrix.sumIntoValues(rowLID,cLIDs,numCols,vals,true,true);
525  }
526  });
527  }
528 
529  }
530  ghostedMatrix->fillComplete(ownedSourceMap,ownedTargetMap);
531  ownedMatrix->resumeFill();
532  ownedMatrix->doExport(*ghostedMatrix,*exporter,Tpetra::ADD);
533  ownedMatrix->fillComplete(ownedSourceMap,ownedTargetMap);
534 
535  return ownedMatrix;
536  }
537 
538 }
539 
540 #endif
Teuchos::RCP< CrsMatrixType > getTpetraMatrix(int i, int j) const
bool is_null(const boost::shared_ptr< T > &p)
virtual Teuchos::RCP< const ExportType > getGhostedExport(int j) const
get exporter for converting an overalapped object to a &quot;normal&quot; object
void setup(const panzer::BasisDescriptor &targetBasis, const panzer::IntegrationDescriptor &integrationDescriptor, const Teuchos::RCP< const Teuchos::MpiComm< int >> &comm, const Teuchos::RCP< const panzer::ConnManager > &connManager, const std::vector< std::string > &elementBlockNames, const Teuchos::RCP< panzer::WorksetContainer > worksetContainer=Teuchos::null)
Setup base objects for L2 Projections - requires target scalar basis and creates worksets if not supp...
const Kokkos::View< const panzer::LocalOrdinal *, Kokkos::LayoutRight, PHX::Device > getElementLIDs(panzer::LocalOrdinal localElmtId) const
virtual int getElementBlockGIDCount(const std::size_t &blockIndex) const =0
How any GIDs are associate with a particular element block.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< panzer::GlobalIndexer > getTargetGlobalIndexer() const
Returns the target global indexer. Will be null if setup() has not been called.
virtual void getOwnedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of indices owned by this processor.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const std::string & getType() const
Get type of basis.
virtual int getFieldNum(const std::string &str) const =0
Get the number used for access to this field.
virtual void getOwnedAndGhostedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of owned and ghosted indices for this processor.
virtual Teuchos::RCP< const MapType > getMap(int i) const
get the map from the matrix
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildRHSMatrix(const panzer::GlobalIndexer &sourceDOFManager, const Teuchos::RCP< const Tpetra::Map< panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device >>> &ownedSourceMap, const std::string &sourceFieldName, const panzer::BasisDescriptor &sourceBasisDescriptor, const int vectorOrGradientDirectionIndex=-1)
Allocates, fills and returns a rectangular matrix for L2 projection of a scalar field, one dimension of gradient (for hgrad basis), or one dimension of a vector field onto the target scalar basis. If you wish to project all values of a vector field or all the gradients of a scalar field, then you will need three separate RHS matrices to form the RHS for each independently. The vectors must be independent Tpetra vectors to solve multiple right hand sides with the linear solver.
Teuchos::RCP< CrsMatrixType > getGhostedTpetraMatrix(int i, int j) const
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 void getElementGIDs(panzer::LocalOrdinal localElmtId, std::vector< panzer::GlobalOrdinal > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
Teuchos::RCP< Tpetra::MultiVector< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildInverseLumpedMassMatrix()
Allocates, fills and returns a Tpetra::MultiVector containing the inverse lumped mass matrix values a...
#define TEUCHOS_ASSERT(assertion_test)
Kokkos::View< const int *, PHX::Device > offsets
Workset size is set to the total number of local elements in the MPI process.
Teuchos::RCP< Tpetra::CrsMatrix< double, panzer::LocalOrdinal, panzer::GlobalOrdinal, Kokkos::Compat::KokkosDeviceWrapperNode< PHX::Device > > > buildMassMatrix(bool use_lumping=false, const std::unordered_map< std::string, double > *elementBlockMultipliers=nullptr)
Allocates, fills and returns a mass matrix for L2 projection onto a target basis. ...