FEI Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
fei_Matrix_Local.hpp
Go to the documentation of this file.
1 /*--------------------------------------------------------------------*/
2 /* Copyright 2007 Sandia Corporation. */
3 /* Under the terms of Contract DE-AC04-94AL85000, there is a */
4 /* non-exclusive license for use of this work by or on behalf */
5 /* of the U.S. Government. Export of this program may require */
6 /* a license from the United States Government. */
7 /*--------------------------------------------------------------------*/
8 
9 #ifndef _fei_Matrix_Local_hpp_
10 #define _fei_Matrix_Local_hpp_
11 
12 #include <fei_SharedPtr.hpp>
13 #include <fei_MatrixGraph.hpp>
14 #include <fei_Matrix.hpp>
15 #include <fei_SparseRowGraph.hpp>
16 
17 #include <vector>
18 
19 namespace fei {
20 class Matrix_Local : public fei::Matrix {
21  public:
24 
25  virtual ~Matrix_Local();
26 
29  bool blockEntry);
30 
31  const char* typeName();
32 
35  int parameters(const fei::ParameterSet& paramset);
36 
39  int parameters(int numParams, const char* const* paramStrings);
40 
43 
46 
49  int getGlobalNumRows() const;
50 
53  int getLocalNumRows() const;
54 
60  int getRowLength(int row, int& length) const;
61 
63  int putScalar(double scalar);
64 
74  int copyOutRow(int row, int len, double* coefs, int* indices) const;
75 
88  int sumIn(int numRows, const int* rows,
89  int numCols, const int* cols,
90  const double* const* values,
91  int format=0);
92 
105  int copyIn(int numRows, const int* rows,
106  int numCols, const int* cols,
107  const double* const* values,
108  int format=0);
109 
125  int sumInFieldData(int fieldID,
126  int idType,
127  int rowID,
128  int colID,
129  const double* const* data,
130  int format=0);
131 
149  int sumInFieldData(int fieldID,
150  int idType,
151  int rowID,
152  int colID,
153  const double* data,
154  int format=0);
155 
165  int sumIn(int blockID, int connectivityID,
166  const double* const* values,
167  int format=0);
168 
173  int globalAssemble();
174 
177  int multiply(fei::Vector* x,
178  fei::Vector* y);
179 
180  void setCommSizes();
181 
187  int gatherFromOverlap(bool accumulate = true);
188 
204  int writeToFile(const char* filename,
205  bool matrixMarketFormat=true);
206 
219  int writeToStream(FEI_OSTREAM& ostrm,
220  bool matrixMarketFormat=true);
221 
224  bool usingBlockEntryStorage();
225 
230  void markState();
231 
236  bool changedSinceMark();
237 
238  const std::vector<int>& getRowNumbers() const;
239 
240  const std::vector<int>& getRowOffsets() const;
241 
242  const std::vector<int>& getColumnIndices() const;
243 
244  const std::vector<double>& getCoefs() const;
245 
246  private:
247  int getRowIndex(int rowNumber) const;
248 
249  int giveToMatrix(int numRows, const int* rows,
250  int numCols, const int* cols,
251  const double* const* values,
252  bool sumInto, int format);
253 
256  std::vector<double> coefs_;
258  std::vector<double> work_data1D_;
259  std::vector<const double*> work_data2D_;
260 };//class Matrix_Local
261 }//namespace fei
262 
263 #endif
264 
Matrix_Local(fei::SharedPtr< fei::MatrixGraph > matrixGraph, fei::SharedPtr< fei::SparseRowGraph > sparseRowGraph)
const std::vector< int > & getRowOffsets() const
int getRowIndex(int rowNumber) const
int sumInFieldData(int fieldID, int idType, int rowID, int colID, const double *const *data, int format=0)
int copyIn(int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)
int sumIn(int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=0)
int getLocalNumRows() const
fei::SharedPtr< fei::SparseRowGraph > sparseRowGraph_
std::vector< double > coefs_
int multiply(fei::Vector *x, fei::Vector *y)
std::vector< const double * > work_data2D_
int giveToMatrix(int numRows, const int *rows, int numCols, const int *cols, const double *const *values, bool sumInto, int format)
int parameters(const fei::ParameterSet &paramset)
const std::vector< int > & getRowNumbers() const
const std::vector< int > & getColumnIndices() const
std::vector< double > work_data1D_
int gatherFromOverlap(bool accumulate=true)
const char * typeName()
std::string filename
#define FEI_OSTREAM
Definition: fei_iosfwd.hpp:24
int writeToStream(FEI_OSTREAM &ostrm, bool matrixMarketFormat=true)
int writeToFile(const char *filename, bool matrixMarketFormat=true)
int getGlobalNumRows() const
int putScalar(double scalar)
const std::vector< double > & getCoefs() const
void setMatrixGraph(fei::SharedPtr< fei::MatrixGraph > matrixGraph)
fei::SharedPtr< fei::MatrixGraph > matrixGraph_
int copyOutRow(int row, int len, double *coefs, int *indices) const
static fei::SharedPtr< fei::Matrix > create_Matrix_Local(fei::SharedPtr< fei::MatrixGraph > matrixGraph, bool blockEntry)
fei::SharedPtr< fei::MatrixGraph > getMatrixGraph() const
int getRowLength(int row, int &length) const