Actual source code: ematrix.c
1: /*
2: Wrappers for esi::petsc::Matrix ESI implementation
3: */
5: #include esi/petsc/matrix.h
7: esi::petsc::Matrix<double,int>::Matrix(esi::IndexSpace<int> *inrmap,esi::IndexSpace<int> *incmap)
8: {
9: int ierr,rn,rN,cn,cN;
10: MPI_Comm *icomm;
12: inrmap->getRunTimeModel("MPI",reinterpret_cast<void *&>(icomm));
13: inrmap->getLocalSize(rn);
14: inrmap->getGlobalSize(rN);
15: incmap->getLocalSize(cn);
16: incmap->getGlobalSize(cN);
17: MatCreate(*icomm,rn,cn,rN,cN,&this->mat);if (ierr) return;
18: PetscObjectSetOptionsPrefix((PetscObject)this->mat,"esi_");
19: MatSetFromOptions(this->mat);
21: this->rmap = inrmap;
22: this->cmap = incmap;
23: inrmap->addReference();
24: incmap->addReference();
26: PetscObjectGetComm((PetscObject)this->mat,&this->comm);if (ierr) return;
27: }
30: esi::petsc::Matrix<double,int>::Matrix(Mat imat)
31: {
32: int m,n,M,N,ierr;
34: this->mat = imat;
35:
36: PetscObjectGetComm((PetscObject)this->mat,&this->comm);if (ierr) return;
37: MatGetLocalSize(mat,&m,&n);if (ierr) return;
38: MatGetSize(mat,&M,&N);if (ierr) return;
39: this->rmap = new esi::petsc::IndexSpace<int>(this->comm,m,M);
40: this->cmap = new esi::petsc::IndexSpace<int>(this->comm,n,N);
41: }
44: esi::petsc::Matrix<double,int>::~Matrix()
45: {
47: if (this->mat) {MatDestroy(this->mat);if (ierr) return;}
48: if (this->rmap) {this->rmap->deleteReference();if (ierr) return;}
49: if (this->cmap) {this->cmap->deleteReference();if (ierr) return;}
50: }
52: esi::ErrorCode esi::petsc::Matrix<double,int>::getInterface(const char* name, void *& iface)
53: {
54: PetscTruth flg;
56: if (!PetscStrcmp(name,"esi::Object",&flg),flg){
57: iface = (void *) (esi::Object *) this;
58: } else if (!PetscStrcmp(name,"esi::Operator",&flg),flg){
59: iface = (void *) (esi::Operator<double,int> *) this;
60: } else if (!PetscStrcmp(name,"esi::MatrixData",&flg),flg){
61: iface = (void *) (esi::MatrixData<int> *) this;
62: } else if (!PetscStrcmp(name,"esi::MatrixRowReadAccess",&flg),flg){
63: iface = (void *) (esi::MatrixRowReadAccess<double,int> *) this;
64: } else if (!PetscStrcmp(name,"esi::MatrixRowWriteAccess",&flg),flg){
65: iface = (void *) (esi::MatrixRowWriteAccess<double,int> *) this;
66: } else if (!PetscStrcmp(name,"esi::petsc::Matrix",&flg),flg){
67: iface = (void *) (esi::petsc::Matrix<double,int> *) this;
68: } else if (!PetscStrcmp(name,"Mat",&flg),flg){
69: iface = (void *) this->mat;
70: } else {
71: iface = 0;
72: }
73: return 0;
74: }
76: esi::ErrorCode esi::petsc::Matrix<double,int>::getInterfacesSupported(esi::Argv * list)
77: {
78: list->appendArg("esi::Object");
79: list->appendArg("esi::Operator");
80: list->appendArg("esi::MatrixData");
81: list->appendArg("esi::MatrixRowReadAccess");
82: list->appendArg("esi::MatrixRowWriteAccess");
83: list->appendArg("esi::petsc::Matrix");
84: list->appendArg("Mat");
85: return 0;
86: }
89: esi::ErrorCode esi::petsc::Matrix<double,int>::apply( esi::Vector<double,int> &xx,esi::Vector<double,int> &yy)
90: {
92: Vec py,px;
94: yy.getInterface("Vec",reinterpret_cast<void*&>(py));
95: xx.getInterface("Vec",reinterpret_cast<void*&>(px));
97: return MatMult(this->mat,px,py);
98: }
100: esi::ErrorCode esi::petsc::Matrix<double,int>::setup()
101: {
103: MatAssemblyBegin(this->mat,MAT_FINAL_ASSEMBLY);
104: return MatAssemblyEnd(this->mat,MAT_FINAL_ASSEMBLY);
105: }
107: esi::ErrorCode esi::petsc::Matrix<double,int>::getIndexSpaces(esi::IndexSpace<int>*& rowIndexSpace, esi::IndexSpace<int>*& colIndexSpace)
108: {
109: rowIndexSpace = this->rmap;
110: colIndexSpace = this->cmap;
111: return 0;
112: }
114: esi::ErrorCode esi::petsc::Matrix<double,int>::isLoaded(bool &State)
115: {
116: return MatAssembled(this->mat,(PetscTruth *)&State);
117: }
119: esi::ErrorCode esi::petsc::Matrix<double,int>::isAllocated(bool &State)
120: {
121: State = 1;
122: return 0;
123: }
125: esi::ErrorCode esi::petsc::Matrix<double,int>::loadComplete()
126: {
127: return this->setup();
128: }
130: esi::ErrorCode esi::petsc::Matrix<double,int>::allocate(int rowlengths[])
131: {
132: return 0;
133: }
135: esi::ErrorCode esi::petsc::Matrix<double,int>::getDiagonal(esi::Vector<double,int> &diagVector)
136: {
138: Vec py;
140: diagVector.getInterface("Vec",reinterpret_cast<void*&>(py));
141: return MatGetDiagonal(this->mat,py);
142: }
144: esi::ErrorCode esi::petsc::Matrix<double,int>::getGlobalSizes(int& rows, int& columns)
145: {
146: return MatGetSize(this->mat,&rows,&columns);
147: }
149: esi::ErrorCode esi::petsc::Matrix<double,int>::getLocalSizes(int& rows, int& columns)
150: {
151: return MatGetLocalSize(this->mat,&rows,&columns);
152: }
154: esi::ErrorCode esi::petsc::Matrix<double,int>::getRowNonzeros(int row,int &length)
155: {
156: int MatGetRow(this->mat,row,&length,PETSC_NULL,PETSC_NULL);
157: return MatRestoreRow(this->mat,row,&length,PETSC_NULL,PETSC_NULL);
158: }
160: esi::ErrorCode esi::petsc::Matrix<double,int>::setRowLength(int row,int length)
161: {
163: return 1;
164: }
166: esi::ErrorCode esi::petsc::Matrix<double,int>::copyOutRow(int row, double* coefs, int* colIndices,int alength,int &length)
167: {
168: int *col;
169: PetscScalar *values;
170:
171: int MatGetRow(this->mat,row,&length,&col,&values);
172: if (length > alength) SETERRQ(1,"Not enough room for values");
173: PetscMemcpy(coefs,values,length*sizeof(PetscScalar));
174: PetscMemcpy(colIndices,col,length*sizeof(int));
175: MatRestoreRow(this->mat,row,&length,&col,&values);
176: return(0);
177: }
179: esi::ErrorCode esi::petsc::Matrix<double,int>::getRow(int row, int& length, double*& coefs, int*& colIndices)
180: {
181: return MatGetRow(this->mat,row,&length,&colIndices,&coefs);
182: }
184: esi::ErrorCode esi::petsc::Matrix<double,int>::getRowCoefs(int row, int& length, double*& coefs)
185: {
186: return MatGetRow(this->mat,row,&length,PETSC_NULL,&coefs);
187: }
189: esi::ErrorCode esi::petsc::Matrix<double,int>::getRowIndices(int row, int& length, int*& colIndices)
190: {
191: return MatGetRow(this->mat,row,&length,&colIndices,PETSC_NULL);
192: }
194: esi::ErrorCode esi::petsc::Matrix<double,int>::restoreRow(int row, int& length, double*& coefs, int*& colIndices)
195: {
196: return MatRestoreRow(this->mat,row,&length,&colIndices,&coefs);
197: }
199: esi::ErrorCode esi::petsc::Matrix<double,int>::restoreRowCoefs(int row, int& length, double*& coefs)
200: {
201: return MatRestoreRow(this->mat,row,&length,PETSC_NULL,&coefs);
202: }
204: esi::ErrorCode esi::petsc::Matrix<double,int>::restoreRowIndices(int row, int& length, int*& colIndices)
205: {
206: return MatRestoreRow(this->mat,row,&length,&colIndices,PETSC_NULL);
207: }
209: esi::ErrorCode esi::petsc::Matrix<double,int>::copyIntoRow(int row,double* coefs, int* colIndices, int length)
210: {
211: return MatSetValues(this->mat,1,&row,length,colIndices,coefs,INSERT_VALUES);
212: }
214: esi::ErrorCode esi::petsc::Matrix<double,int>::sumIntoRow(int row, double* coefs, int* colIndices,int length)
215: {
216: return MatSetValues(this->mat,1,&row,length,colIndices,coefs,ADD_VALUES);
217: }
219: esi::ErrorCode esi::petsc::Matrix<double,int>::rowMax(int row, double &result)
220: {
221: int ierr,length,i;
222: PetscScalar *values;
224: MatGetRow(this->mat,row,&length,PETSC_NULL,&values);
225: if (values) {
226: result = values[0];
227: for (i=1; i<length; i++) result = PetscMax(result,values[i]);
228: }
229: MatRestoreRow(this->mat,row,&length,PETSC_NULL,&values);
230: return(0);
231: }
233: esi::ErrorCode esi::petsc::Matrix<double,int>::rowMin(int row, double &result)
234: {
235: int ierr,length,i;
236: PetscScalar *values;
238: MatGetRow(this->mat,row,&length,PETSC_NULL,&values);
239: if (values) {
240: result = values[0];
241: for (i=1; i<length; i++) result = PetscMin(result,values[i]);
242: }
243: MatRestoreRow(this->mat,row,&length,PETSC_NULL,&values);
244: return(0);
245: }
247: esi::ErrorCode esi::petsc::Matrix<double,int>::getRowSum(esi::Vector<double,int>& rowSumVector)
248: {
249: int i,ierr,rstart,rend,length,j;
250: PetscScalar *values,sum;
251: Vec py;
253: rowSumVector.getInterface("Vec",reinterpret_cast<void*&>(py));
255: MatGetOwnershipRange(this->mat,&rstart,&rend);
256: for ( i=rstart; i<rend; i++) {
257: MatGetRow(this->mat,i,&length,PETSC_NULL,&values);
258: sum = 0.0;
259: for ( j=0; j<length; j++ ) {
260: sum += values[j];
261: }
262: MatRestoreRow(this->mat,i,&length,PETSC_NULL,&values);
263: VecSetValues(py,1,&i,&sum,INSERT_VALUES);
264: }
265: VecAssemblyBegin(py);
266: VecAssemblyEnd(py);
268: return 0;
269: }
271: /* ------------------------------------------------------------------------------------------------*/
273: ::esi::ErrorCode esi::petsc::Matrix<double,int>::Factory::create(::esi::IndexSpace<int>&irmap,::esi::IndexSpace<int>&icmap,::esi::Operator<double,int>*&v)
274: {
275: v = new esi::petsc::Matrix<double,int>(&irmap,&icmap);
276: return 0;
277: };
279: /* ::esi::petsc::OperatorFactory<double,int> OFInstForIntel64CompilerBug; */
281: EXTERN_C_BEGIN
282: ::esi::Operator<double,int>::Factory *create_esi_petsc_operatorfactory(void)
283: {
284: return dynamic_cast< ::esi::Operator<double,int>::Factory *>(new esi::petsc::Matrix<double,int>::Factory);
285: }
286: EXTERN_C_END
289: #if defined(PETSC_HAVE_TRILINOS)
290: #include esi/petsc/matrix.h
291: #define PETRA_MPI /* used by Ptera to indicate MPI code */
292: #include "Petra_ESI_CRS_Matrix.h"
294: /*
295: This class is the same as the Petra_ESI_CRS_Matrix class except it puts values into the Petra_CRS_Grap()
296: */
297: template<class Scalar,class Ordinal> class MyPetra_ESI_CRS_Matrix : public virtual Petra_ESI_CRS_Matrix<Scalar,Ordinal>
298: {
299: public:
301: MyPetra_ESI_CRS_Matrix(Petra_DataAccess CV,const Petra_CRS_Graph& graph) : Petra_ESI_Object(), Petra_RDP_CRS_Matrix(CV, graph), Petra_ESI_CRS_Matrix<Scalar,Ordinal>(CV, graph){graph_ = (Petra_CRS_Graph*)&graph;SetStaticGraph(false);};
303: virtual ~MyPetra_ESI_CRS_Matrix(void) { };
305: virtual esi::ErrorCode copyIntoRow(Ordinal row, Scalar* coefs, Ordinal* colIndices, Ordinal length)
307: graph_->InsertGlobalIndices(row, length, colIndices);CHKERRQ(((ierr == 1) ? 0 : ierr));
308: Petra_ESI_CRS_Matrix<Scalar,Ordinal>::copyIntoRow(row,coefs,colIndices,length);
309: // this->setup();
310: return 0;
311: }
313: Petra_CRS_Graph *graph_;
314: };
317: template<class Scalar,class Ordinal> class Petra_ESI_CRS_OperatorFactory : public virtual ::esi::OperatorFactory<Scalar,Ordinal>
318: {
319: public:
321: // constructor
322: Petra_ESI_CRS_OperatorFactory(void){};
323:
324: // Destructor.
325: virtual ~Petra_ESI_CRS_OperatorFactory(void){};
328: // Construct a Operator
329: virtual ::esi::ErrorCode getOperator(::esi::IndexSpace<Ordinal>&rmap,::esi::IndexSpace<Ordinal>&cmap,::esi::Operator<Scalar,Ordinal>*&v)
330: {
331: int ierr;
332: Petra_Map *rowmap,*colmap;
333: rmap.getInterface("Petra_Map",reinterpret_cast<void *&>(rowmap));
334: if (!rowmap) SETERRQ(1,"Petra requires all IndexSpaces be Petra_ESI_IndexSpaces");
335: cmap.getInterface("Petra_Map",reinterpret_cast<void *&>(colmap));
336: if (!colmap) SETERRQ(1,"Petra requires all IndexSpaces be Petra_ESI_IndexSpaces");
337: Petra_CRS_Graph *graph = new Petra_CRS_Graph(Copy,*(Petra_BlockMap*)rowmap,*(Petra_BlockMap*)colmap,200);
338: rmap.addReference();
339: cmap.addReference();
340: v = new MyPetra_ESI_CRS_Matrix<Scalar,Ordinal>(Copy,*graph);
341: return 0;
342: };
343: };
345: EXTERN_C_BEGIN
346: ::esi::OperatorFactory<double,int> *create_petra_esi_operatorfactory(void)
347: {
348: return dynamic_cast< ::esi::OperatorFactory<double,int> *>(new Petra_ESI_CRS_OperatorFactory<double,int>);
349: }
350: EXTERN_C_END
351: #endif
355: