Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_AdaptivityUtils.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
13 
14 #include "Epetra_IntVector.h"
15 #include "Epetra_Import.h"
16 
18  const Epetra_Comm & Comm,
19  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_row_basis,
20  std::vector<int> & myRowGidOffsets)
21 {
22  // build row offsets
23  int totalSz = 0;
24  std::vector<int> myRowOffsets(per_dof_row_basis.size());
25  for(std::size_t i=0;i<per_dof_row_basis.size();i++) {
26  myRowOffsets[i] = totalSz;
27  totalSz += per_dof_row_basis[i]->size();
28  }
29 
30  // build row map
31  Teuchos::RCP<Epetra_Map> rowMap = Teuchos::rcp(new Epetra_Map(-1,totalSz,0,Comm));
32 
33  myRowGidOffsets.resize(per_dof_row_basis.size());
34  for(std::size_t i=0;i<myRowGidOffsets.size();i++)
35  myRowGidOffsets[i] = rowMap->GID(myRowOffsets[i]);
36 
37  return rowMap;
38 }
39 
41  const Epetra_Comm & Comm,
42  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_row_basis)
43 {
44  // build row offsets
45  int totalSz = 0;
46  for(std::size_t i=0;i<per_dof_row_basis.size();i++)
47  totalSz += per_dof_row_basis[i]->size();
48 
49  // build row map
50  return Teuchos::rcp(new Epetra_Map(-1,totalSz,0,Comm));
51 }
52 
54  const Epetra_CrsGraph & determGraph,
55  const std::vector<int> & myRowGidOffsets,
56  std::vector<int> & myColGidOffsets)
57 {
58  // build global distributed offsets index
59  Epetra_IntVector adaptRowGids(determGraph.RowMap()); // uniquely defined row GIDs that will
60  // be used to determine col GIDs
61  for(std::size_t i=0;i<myRowGidOffsets.size();i++)
62  adaptRowGids[i] = myRowGidOffsets[i];
63 
64  // perform global communication to determine offsets
65  Epetra_Import importer(determGraph.ColMap(),determGraph.RowMap());
66  Epetra_IntVector adaptColGids(determGraph.ColMap());
67  adaptColGids.Import(adaptRowGids,importer,Insert);
68 
69  // copy offsets into std::vector
70  myColGidOffsets.resize(adaptColGids.MyLength());
71  for(std::size_t i=0;i<myColGidOffsets.size();i++)
72  myColGidOffsets[i] = adaptColGids[i];
73 }
74 
76  const Epetra_CrsGraph & determGraph,
77  const Teuchos::RCP<const Stokhos::ProductBasis<int,double> > & masterBasis,
78  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_row_basis,
79  std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_col_basis)
80 {
81  using Teuchos::RCP;
82 
83  const Epetra_BlockMap & rowMap = determGraph.RowMap();
84  const Epetra_BlockMap & colMap = determGraph.ColMap();
85 
86  // this is block size
87  int numStochDim = masterBasis->dimension();
88 
89  // build blocked maps
90  Epetra_BlockMap blkRowMap(-1,rowMap.NumMyElements(),rowMap.MyGlobalElements(),numStochDim,0,rowMap.Comm());
91  Epetra_BlockMap blkColMap(-1,colMap.NumMyElements(),colMap.MyGlobalElements(),numStochDim,0,colMap.Comm());
92  // construct int vectors to determine order
93  Epetra_IntVector stochRowOrders(blkRowMap),
94  stochColOrders(blkColMap); // to build built by Import
95 
96  // loop over row degrees of freedom building Row Orders vector
97  for(std::size_t dof=0;dof<per_dof_row_basis.size();dof++) {
98  RCP<const Stokhos::ProductBasis<int,double> > rowBasis = per_dof_row_basis[dof];
99  TEUCHOS_TEST_FOR_EXCEPTION(rowBasis->dimension()!=masterBasis->dimension(),std::invalid_argument,
100  "Stokhos::adapt_utils::buildColBasisFunctions: Row basis must match dimension of master basis!");
101 
103  = rowBasis->getCoordinateBases();
104 
105  TEUCHOS_TEST_FOR_EXCEPTION(onedBasis.size()!=numStochDim,std::logic_error,
106  "Stokhos::adapt_utils::buildColBasisFunctions: Wrong number of dimensions from row basis!");
107 
108  // fill stochastic orders vector
109  for(int i=0;i<numStochDim;i++)
110  stochRowOrders[i+dof*numStochDim] = onedBasis[i]->order();
111  }
112 
113  // do communication to determine row maps
114  Epetra_Import importer(blkColMap,blkRowMap);
115  stochColOrders.Import(stochRowOrders,importer,Insert);
116 
118  = masterBasis->getCoordinateBases();
119 
120  // now construct column basis functions
121  std::vector<int> polyOrder(numStochDim,0);
122  per_dof_col_basis.resize(blkColMap.NumMyElements());
123  for(std::size_t col=0;col<per_dof_col_basis.size();col++) {
124  int colOffset = numStochDim*col;
125 
126  // extract polynomial order for this column
127  for(int o=0;o<numStochDim;o++)
128  polyOrder[o] = stochColOrders[colOffset+o];
129 
131  for(Teuchos::Ordinal dim=0;dim<oneDMasterBasis.size();dim++) {
132  RCP<const OneDOrthogPolyBasis<int,double> > oneDBasis = oneDMasterBasis[dim];
133 
134  newBasisArray[dim] = oneDBasis->cloneWithOrder(polyOrder[dim]);
135  }
136 
137  per_dof_col_basis[col] = Teuchos::rcp(
139  }
140 }
141 
143  const Epetra_CrsGraph & determGraph,
144  const Teuchos::RCP<const Stokhos::ProductBasis<int,double> > & masterBasis,
145  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_row_basis,
146  bool onlyUseLinear,
147  int kExpOrder)
148 {
149  std::vector<int> myRowGidOffsets, myColGidOffsets;
150 
151  return buildAdaptedGraph(determGraph,masterBasis,per_dof_row_basis,myRowGidOffsets,myColGidOffsets,onlyUseLinear,kExpOrder);
152 }
153 
155  const Epetra_CrsGraph & determGraph,
156  const Teuchos::RCP<const Stokhos::ProductBasis<int,double> > & masterBasis,
157  const std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > & per_dof_row_basis,
158  std::vector<int> & myRowGidOffsets,std::vector<int> & myColGidOffsets,
159  bool onlyUseLinear,
160  int kExpOrder)
161 {
162  TEUCHOS_TEST_FOR_EXCEPTION(int(per_dof_row_basis.size())!=determGraph.NumMyRows(),std::logic_error,
163  "Stokhos::adapt_utils::buildAdaptedGraph: per_dof_row_basis.size()!=determGraph.NumMyRows()");
164 
165  myRowGidOffsets.clear();
166  myColGidOffsets.clear();
167 
168  std::vector<Teuchos::RCP<const Stokhos::ProductBasis<int,double> > > per_dof_col_basis;
170 
171  // build basis functions associated with the columns
172  buildColBasisFunctions(determGraph,masterBasis,per_dof_row_basis,per_dof_col_basis);
173 
174  // build row map, and row and column offsets.
175  rowMap = buildAdaptedRowMapAndOffsets(determGraph.Comm(),per_dof_row_basis,myRowGidOffsets);
176  buildAdaptedColOffsets(determGraph,myRowGidOffsets,myColGidOffsets);
177 
179 
181  if(kExpOrder<0)
182  Cijk = masterBasis->computeTripleProductTensor();
183  else
184  Cijk = masterBasis->computeLinearTripleProductTensor();
185 
186  // iterate over nonzero structure of graph
187  int maxNNZ = determGraph.MaxNumNonzeros();
188  std::vector<int> determGraphCols(maxNNZ);
189  std::vector<int> graphCols;
190  for(int lRID=0;lRID<determGraph.NumMyRows();lRID++) {
191  int gRID = determGraph.GRID(lRID);
192  int numIndices = -1;
193  int rowOffsetIndex = myRowGidOffsets[lRID];
194 
195  // grab row nonzero structure
196  determGraph.ExtractGlobalRowCopy(gRID,maxNNZ,numIndices,&determGraphCols[0]);
197 
198  for(int i=0;i<numIndices;i++) {
199  int gCID = determGraphCols[i];
200  int lCID = determGraph.LCID(gCID);
201  int colOffsetIndex = myColGidOffsets[lCID];
202 
203  Stokhos::BasisInteractionGraph interactGraph(*masterBasis,*Cijk,*per_dof_row_basis[lRID],
204  *per_dof_col_basis[lCID],
205  onlyUseLinear,kExpOrder);
206  for(std::size_t basisRow=0;basisRow<interactGraph.rowCount();basisRow++) {
207  const std::vector<std::size_t> & basisCols = interactGraph.activeIndices(basisRow);
208  graphCols.resize(basisCols.size());
209  for(std::size_t g=0;g<basisCols.size();g++)
210  graphCols[g] = basisCols[g] + colOffsetIndex;
211 
212  graph->InsertGlobalIndices(basisRow+rowOffsetIndex,graphCols.size(),&graphCols[0]);
213  }
214  }
215  }
216 
217  graph->FillComplete();
218 
219  return graph;
220 }
Teuchos::RCP< Epetra_Map > buildAdaptedRowMap(const Epetra_Comm &Comm, const std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > &per_dof_row_basis)
const Epetra_Comm & Comm() const
int MyGlobalElements(int *MyGlobalElementList) const
void buildColBasisFunctions(const Epetra_CrsGraph &determGraph, const Teuchos::RCP< const Stokhos::ProductBasis< int, double > > &masterBasis, const std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > &per_dof_row_basis, std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > &per_dof_col_basis)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
int NumMyRows() const
const Epetra_BlockMap & ColMap() const
const std::vector< std::size_t > & activeIndices(std::size_t i) const
Grab active indicies in graph for row i.
int InsertGlobalIndices(int_type GlobalRow, int NumIndices, int_type *Indices)
int MaxNumNonzeros() const
int NumMyElements() const
Teuchos::RCP< Epetra_CrsGraph > buildAdaptedGraph(const Epetra_CrsGraph &determGraph, const Teuchos::RCP< const Stokhos::ProductBasis< int, double > > &masterBasis, const std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > &per_dof_row_basis, bool onlyUseLinear=false, int kExpOrder=-1)
int LCID(int GCID_in) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int GID(int LID) const
Teuchos::RCP< Epetra_Map > buildAdaptedRowMapAndOffsets(const Epetra_Comm &Comm, const std::vector< Teuchos::RCP< const Stokhos::ProductBasis< int, double > > > &per_dof_row_basis, std::vector< int > &myRowGidOffsets)
const Epetra_BlockMap & RowMap() const
void buildAdaptedColOffsets(const Epetra_CrsGraph &determGraph, const std::vector< int > &myRowGidOffsets, std::vector< int > &myColGidOffsets)
const Epetra_Comm & Comm() const
int GRID(int LRID_in) const
std::size_t rowCount() const
What is the number of rows.
size_type size() const
int ExtractGlobalRowCopy(int_type Row, int LenOfIndices, int &NumIndices, int_type *Indices) const
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)