xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 17699dbb92d734bf05dd659d53b72cde1e635442)
1cb512458SBarry Smith #ifndef lint
2*17699dbbSLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.66 1995/10/18 03:18:20 curfman Exp curfman $";
3cb512458SBarry Smith #endif
4289bc588SBarry Smith 
5*17699dbbSLois Curfman McInnes #include "dense.h"
6b16a3bb1SBarry Smith #include "pinclude/plapack.h"
7b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
8289bc588SBarry Smith 
91987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
101987afe7SBarry Smith {
111987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
121987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
131987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
141987afe7SBarry Smith   PLogFlops(2*N-1);
151987afe7SBarry Smith   return 0;
161987afe7SBarry Smith }
171987afe7SBarry Smith 
18c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
19289bc588SBarry Smith {
20c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
21289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
22289bc588SBarry Smith   Scalar       *v = mat->v;
23289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
24c0bbcb79SLois Curfman McInnes   *nz = count; *nzalloc = N; *mem = (int)A->mem;
25fa9ec3f1SLois Curfman McInnes   return 0;
26289bc588SBarry Smith }
27289bc588SBarry Smith 
28289bc588SBarry Smith /* ---------------------------------------------------------------*/
29289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
30289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
31c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
32289bc588SBarry Smith {
33c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
346abc6512SBarry Smith   int          info;
35289bc588SBarry Smith   if (!mat->pivots) {
3648d91487SBarry Smith     mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
37c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
38289bc588SBarry Smith   }
39289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
40ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
41c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
42289bc588SBarry Smith   return 0;
43289bc588SBarry Smith }
44c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
45289bc588SBarry Smith {
46289bc588SBarry Smith   int ierr;
47c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
48289bc588SBarry Smith   return 0;
49289bc588SBarry Smith }
50c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
51289bc588SBarry Smith {
5249d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
53289bc588SBarry Smith }
54c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
55289bc588SBarry Smith {
56289bc588SBarry Smith   int ierr;
57c0bbcb79SLois Curfman McInnes   ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr);
58289bc588SBarry Smith   return 0;
59289bc588SBarry Smith }
60c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
61289bc588SBarry Smith {
6249d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
63289bc588SBarry Smith }
64c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
65289bc588SBarry Smith {
66c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
676abc6512SBarry Smith   int           info;
68752f5567SLois Curfman McInnes   if (mat->pivots) {
69752f5567SLois Curfman McInnes     PETSCFREE(mat->pivots);
70c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
71752f5567SLois Curfman McInnes     mat->pivots = 0;
72752f5567SLois Curfman McInnes   }
73289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
74ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
75c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
76289bc588SBarry Smith   return 0;
77289bc588SBarry Smith }
78289bc588SBarry Smith 
79c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
80289bc588SBarry Smith {
81c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
826abc6512SBarry Smith   int          one = 1, info;
836abc6512SBarry Smith   Scalar       *x, *y;
84289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
85416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
86c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
8748d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
88289bc588SBarry Smith   }
89c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
9048d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
91289bc588SBarry Smith   }
92ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
93ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
94289bc588SBarry Smith   return 0;
95289bc588SBarry Smith }
96c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
97da3a660dSBarry Smith {
98c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
996abc6512SBarry Smith   int          one = 1, info;
1006abc6512SBarry Smith   Scalar       *x, *y;
101da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
102416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
103752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
104da3a660dSBarry Smith   if (mat->pivots) {
10548d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
106da3a660dSBarry Smith   }
107da3a660dSBarry Smith   else {
10848d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
109da3a660dSBarry Smith   }
110ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
111da3a660dSBarry Smith   return 0;
112da3a660dSBarry Smith }
113c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
114da3a660dSBarry Smith {
115c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1166abc6512SBarry Smith   int          one = 1, info,ierr;
1176abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
118da3a660dSBarry Smith   Vec          tmp = 0;
119da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
120da3a660dSBarry Smith   if (yy == zz) {
12178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
122c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
12378b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
124da3a660dSBarry Smith   }
125416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
126752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
127da3a660dSBarry Smith   if (mat->pivots) {
12848d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
129da3a660dSBarry Smith   }
130da3a660dSBarry Smith   else {
13148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
132da3a660dSBarry Smith   }
133ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
134da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
135da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
136da3a660dSBarry Smith   return 0;
137da3a660dSBarry Smith }
138c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
139da3a660dSBarry Smith {
140c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1416abc6512SBarry Smith   int           one = 1, info,ierr;
1426abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
143da3a660dSBarry Smith   Vec           tmp;
144da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
145da3a660dSBarry Smith   if (yy == zz) {
14678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
147c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
14878b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
149da3a660dSBarry Smith   }
150416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
151752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
152da3a660dSBarry Smith   if (mat->pivots) {
15348d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
154da3a660dSBarry Smith   }
155da3a660dSBarry Smith   else {
15648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
157da3a660dSBarry Smith   }
158ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
159da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
160da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
161da3a660dSBarry Smith   return 0;
162da3a660dSBarry Smith }
163289bc588SBarry Smith /* ------------------------------------------------------------------*/
164c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
16520e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
166289bc588SBarry Smith {
167c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
168289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1696abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
170289bc588SBarry Smith 
171289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
172289bc588SBarry Smith     /* this is a hack fix, should have another version without
173289bc588SBarry Smith        the second BLdot */
174bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
175289bc588SBarry Smith   }
176289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
177289bc588SBarry Smith   while (its--) {
178289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
179289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
180f1747703SBarry Smith #if defined(PETSC_COMPLEX)
181f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
182f1747703SBarry Smith            not happy about returning a double complex */
183f1747703SBarry Smith         int    _i;
184f1747703SBarry Smith         Scalar sum = b[i];
185f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
186f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
187f1747703SBarry Smith         }
188f1747703SBarry Smith         xt = sum;
189f1747703SBarry Smith #else
190289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
191f1747703SBarry Smith #endif
192d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
193289bc588SBarry Smith       }
194289bc588SBarry Smith     }
195289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
196289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
197f1747703SBarry Smith #if defined(PETSC_COMPLEX)
198f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
199f1747703SBarry Smith            not happy about returning a double complex */
200f1747703SBarry Smith         int    _i;
201f1747703SBarry Smith         Scalar sum = b[i];
202f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
203f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
204f1747703SBarry Smith         }
205f1747703SBarry Smith         xt = sum;
206f1747703SBarry Smith #else
207289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
208f1747703SBarry Smith #endif
209d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
210289bc588SBarry Smith       }
211289bc588SBarry Smith     }
212289bc588SBarry Smith   }
213289bc588SBarry Smith   return 0;
214289bc588SBarry Smith }
215289bc588SBarry Smith 
216289bc588SBarry Smith /* -----------------------------------------------------------------*/
217c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
218289bc588SBarry Smith {
219c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
220289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
221289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
222289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
22348d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
224289bc588SBarry Smith   return 0;
225289bc588SBarry Smith }
226c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
227289bc588SBarry Smith {
228c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
229289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
230289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
231289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
23248d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
233289bc588SBarry Smith   return 0;
234289bc588SBarry Smith }
235c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
236289bc588SBarry Smith {
237c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
238289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2396abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
240289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
241416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
24248d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
243289bc588SBarry Smith   return 0;
244289bc588SBarry Smith }
245c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
246289bc588SBarry Smith {
247c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
248289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
249289bc588SBarry Smith   int          _One=1;
2506abc6512SBarry Smith   Scalar       _DOne=1.0;
251289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
252289bc588SBarry Smith   VecGetArray(zz,&z);
253416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25448d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
255289bc588SBarry Smith   return 0;
256289bc588SBarry Smith }
257289bc588SBarry Smith 
258289bc588SBarry Smith /* -----------------------------------------------------------------*/
259c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
260289bc588SBarry Smith {
261c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
262289bc588SBarry Smith   Scalar       *v;
263289bc588SBarry Smith   int          i;
264289bc588SBarry Smith   *ncols = mat->n;
265289bc588SBarry Smith   if (cols) {
26678b31e54SBarry Smith     *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols);
26773c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
268289bc588SBarry Smith   }
269289bc588SBarry Smith   if (vals) {
27078b31e54SBarry Smith     *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
271289bc588SBarry Smith     v = mat->v + row;
27273c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
273289bc588SBarry Smith   }
274289bc588SBarry Smith   return 0;
275289bc588SBarry Smith }
276c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
277289bc588SBarry Smith {
27878b31e54SBarry Smith   if (cols) { PETSCFREE(*cols); }
27978b31e54SBarry Smith   if (vals) { PETSCFREE(*vals); }
280289bc588SBarry Smith   return 0;
281289bc588SBarry Smith }
282289bc588SBarry Smith /* ----------------------------------------------------------------*/
283c0bbcb79SLois Curfman McInnes static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n,
284e8d4e0b9SBarry Smith                         int *indexn,Scalar *v,InsertMode addv)
285289bc588SBarry Smith {
286c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
287289bc588SBarry Smith   int          i,j;
288d6dfbf8fSBarry Smith 
289289bc588SBarry Smith   if (!mat->roworiented) {
290dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
291289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
292289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
293289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
294289bc588SBarry Smith         }
295289bc588SBarry Smith       }
296289bc588SBarry Smith     }
297289bc588SBarry Smith     else {
298289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
299289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
300289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
301289bc588SBarry Smith         }
302289bc588SBarry Smith       }
303289bc588SBarry Smith     }
304e8d4e0b9SBarry Smith   }
305e8d4e0b9SBarry Smith   else {
306dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
307e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
308e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
309e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
310e8d4e0b9SBarry Smith         }
311e8d4e0b9SBarry Smith       }
312e8d4e0b9SBarry Smith     }
313289bc588SBarry Smith     else {
314289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
315289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
316289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
317289bc588SBarry Smith         }
318289bc588SBarry Smith       }
319289bc588SBarry Smith     }
320e8d4e0b9SBarry Smith   }
321289bc588SBarry Smith   return 0;
322289bc588SBarry Smith }
323e8d4e0b9SBarry Smith 
324289bc588SBarry Smith /* -----------------------------------------------------------------*/
325c0bbcb79SLois Curfman McInnes static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat)
326289bc588SBarry Smith {
327c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
328289bc588SBarry Smith   int          ierr;
329289bc588SBarry Smith   Mat          newi;
33048d91487SBarry Smith 
331c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi); CHKERRQ(ierr);
332ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
333416022c9SBarry Smith   PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
334289bc588SBarry Smith   *newmat = newi;
335289bc588SBarry Smith   return 0;
336289bc588SBarry Smith }
337289bc588SBarry Smith 
338932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
339932b0c3eSLois Curfman McInnes #include "sysio.h"
340932b0c3eSLois Curfman McInnes 
341932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
342289bc588SBarry Smith {
343932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
344932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
345d636dbe3SBarry Smith   FILE         *fd;
346932b0c3eSLois Curfman McInnes   char         *outputname;
347932b0c3eSLois Curfman McInnes   Scalar       *v;
348932b0c3eSLois Curfman McInnes 
349932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
350932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
351932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetFormat_Private(viewer,&format);
352f72585f2SLois Curfman McInnes   if (format == FILE_FORMAT_INFO) {
353932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
354f72585f2SLois Curfman McInnes   }
355f72585f2SLois Curfman McInnes   else {
356932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
357932b0c3eSLois Curfman McInnes       v = a->v + i;
358932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
359289bc588SBarry Smith #if defined(PETSC_COMPLEX)
360932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
361289bc588SBarry Smith #else
362932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
363289bc588SBarry Smith #endif
364289bc588SBarry Smith       }
3658e37dc09SBarry Smith       fprintf(fd,"\n");
366289bc588SBarry Smith     }
367da3a660dSBarry Smith   }
368932b0c3eSLois Curfman McInnes   fflush(fd);
369289bc588SBarry Smith   return 0;
370289bc588SBarry Smith }
371289bc588SBarry Smith 
372932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
373932b0c3eSLois Curfman McInnes {
374932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
375932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
376932b0c3eSLois Curfman McInnes   Scalar       *v, *anonz;
377932b0c3eSLois Curfman McInnes 
378932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
379932b0c3eSLois Curfman McInnes   col_lens = (int *) PETSCMALLOC( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
380932b0c3eSLois Curfman McInnes   col_lens[0] = MAT_COOKIE;
381932b0c3eSLois Curfman McInnes   col_lens[1] = m;
382932b0c3eSLois Curfman McInnes   col_lens[2] = n;
383932b0c3eSLois Curfman McInnes   col_lens[3] = nz;
384932b0c3eSLois Curfman McInnes 
385932b0c3eSLois Curfman McInnes   /* store lengths of each row and write (including header) to file */
386932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) col_lens[4+i] = n;
387932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr);
388932b0c3eSLois Curfman McInnes 
389932b0c3eSLois Curfman McInnes   /* Possibly should write in smaller increments, not whole matrix at once? */
390932b0c3eSLois Curfman McInnes  /* store column indices (zero start index) */
391932b0c3eSLois Curfman McInnes   ict = 0;
392932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
393932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) col_lens[ict++] = j;
394932b0c3eSLois Curfman McInnes   }
395932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr);
396932b0c3eSLois Curfman McInnes   PETSCFREE(col_lens);
397932b0c3eSLois Curfman McInnes 
398932b0c3eSLois Curfman McInnes   /* store nonzero values */
399932b0c3eSLois Curfman McInnes   anonz = (Scalar *) PETSCMALLOC((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
400932b0c3eSLois Curfman McInnes   ict = 0;
401932b0c3eSLois Curfman McInnes   for ( i=0; i<m; i++ ) {
402932b0c3eSLois Curfman McInnes     v = a->v + i;
403932b0c3eSLois Curfman McInnes     for ( j=0; j<n; j++ ) {
404932b0c3eSLois Curfman McInnes       anonz[ict++] = *v; v += a->m;
405932b0c3eSLois Curfman McInnes     }
406932b0c3eSLois Curfman McInnes   }
407932b0c3eSLois Curfman McInnes   ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr);
408932b0c3eSLois Curfman McInnes   PETSCFREE(anonz);
409932b0c3eSLois Curfman McInnes   return 0;
410932b0c3eSLois Curfman McInnes }
411932b0c3eSLois Curfman McInnes 
412932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
413932b0c3eSLois Curfman McInnes {
414932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
415932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
416932b0c3eSLois Curfman McInnes   PetscObject  vobj = (PetscObject) viewer;
417932b0c3eSLois Curfman McInnes 
418932b0c3eSLois Curfman McInnes   if (!viewer) {
419932b0c3eSLois Curfman McInnes     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
420932b0c3eSLois Curfman McInnes   }
421932b0c3eSLois Curfman McInnes   if (vobj->cookie == VIEWER_COOKIE) {
422932b0c3eSLois Curfman McInnes     if (vobj->type == MATLAB_VIEWER) {
423932b0c3eSLois Curfman McInnes       return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
424932b0c3eSLois Curfman McInnes     }
425932b0c3eSLois Curfman McInnes     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) {
426932b0c3eSLois Curfman McInnes       return MatView_SeqDense_ASCII(A,viewer);
427932b0c3eSLois Curfman McInnes     }
428932b0c3eSLois Curfman McInnes     else if (vobj->type == BINARY_FILE_VIEWER) {
429932b0c3eSLois Curfman McInnes       return MatView_SeqDense_Binary(A,viewer);
430932b0c3eSLois Curfman McInnes     }
431932b0c3eSLois Curfman McInnes   }
432932b0c3eSLois Curfman McInnes   return 0;
433932b0c3eSLois Curfman McInnes }
434289bc588SBarry Smith 
435ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
436289bc588SBarry Smith {
437289bc588SBarry Smith   Mat          mat = (Mat) obj;
438ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
439a5a9c739SBarry Smith #if defined(PETSC_LOG)
440a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
441a5a9c739SBarry Smith #endif
44278b31e54SBarry Smith   if (l->pivots) PETSCFREE(l->pivots);
44378b31e54SBarry Smith   PETSCFREE(l);
444a5a9c739SBarry Smith   PLogObjectDestroy(mat);
4459e25ed09SBarry Smith   PETSCHEADERDESTROY(mat);
446289bc588SBarry Smith   return 0;
447289bc588SBarry Smith }
448289bc588SBarry Smith 
449c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
450289bc588SBarry Smith {
451c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
452d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
453d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
45448b35521SBarry Smith 
455d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
45648b35521SBarry Smith   if (!matout) { /* in place transpose */
45748d91487SBarry Smith     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place");
458d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
459289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
460d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
461d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
462d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
463289bc588SBarry Smith       }
464289bc588SBarry Smith     }
46548b35521SBarry Smith   }
466d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
467d3e5ee88SLois Curfman McInnes     int          ierr;
468d3e5ee88SLois Curfman McInnes     Mat          tmat;
469ec8511deSBarry Smith     Mat_SeqDense *tmatd;
470d3e5ee88SLois Curfman McInnes     Scalar       *v2;
471c0bbcb79SLois Curfman McInnes     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr);
472ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
4730de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
474d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
4750de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
476d3e5ee88SLois Curfman McInnes     }
477d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
478d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
479d3e5ee88SLois Curfman McInnes     *matout = tmat;
48048b35521SBarry Smith   }
481289bc588SBarry Smith   return 0;
482289bc588SBarry Smith }
483289bc588SBarry Smith 
484c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2)
485289bc588SBarry Smith {
486c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
487c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
488289bc588SBarry Smith   int          i;
489289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
490289bc588SBarry Smith   if (mat1->m != mat2->m) return 0;
491289bc588SBarry Smith   if (mat1->n != mat2->n) return 0;
492289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
493289bc588SBarry Smith     if (*v1 != *v2) return 0;
494289bc588SBarry Smith     v1++; v2++;
495289bc588SBarry Smith   }
496289bc588SBarry Smith   return 1;
497289bc588SBarry Smith }
498289bc588SBarry Smith 
499c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
500289bc588SBarry Smith {
501c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
5026abc6512SBarry Smith   int          i, n;
5036abc6512SBarry Smith   Scalar       *x;
504289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
505ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
506289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
507289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
508289bc588SBarry Smith   }
509289bc588SBarry Smith   return 0;
510289bc588SBarry Smith }
511289bc588SBarry Smith 
512c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr)
513289bc588SBarry Smith {
514c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
515da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
516da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
51728988994SBarry Smith   if (ll) {
518da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
519ec8511deSBarry Smith     if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size");
520da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
521da3a660dSBarry Smith       x = l[i];
522da3a660dSBarry Smith       v = mat->v + i;
523da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
524da3a660dSBarry Smith     }
525da3a660dSBarry Smith   }
52628988994SBarry Smith   if (rr) {
527da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
528ec8511deSBarry Smith     if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size");
529da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
530da3a660dSBarry Smith       x = r[i];
531da3a660dSBarry Smith       v = mat->v + i*m;
532da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
533da3a660dSBarry Smith     }
534da3a660dSBarry Smith   }
535289bc588SBarry Smith   return 0;
536289bc588SBarry Smith }
537289bc588SBarry Smith 
538c0bbcb79SLois Curfman McInnes static int MatNorm_SeqDense(Mat A,MatNormType type,double *norm)
539289bc588SBarry Smith {
540c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
541289bc588SBarry Smith   Scalar       *v = mat->v;
542289bc588SBarry Smith   double       sum = 0.0;
543289bc588SBarry Smith   int          i, j;
544289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
545289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
546289bc588SBarry Smith #if defined(PETSC_COMPLEX)
547289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
548289bc588SBarry Smith #else
549289bc588SBarry Smith       sum += (*v)*(*v); v++;
550289bc588SBarry Smith #endif
551289bc588SBarry Smith     }
552289bc588SBarry Smith     *norm = sqrt(sum);
553289bc588SBarry Smith   }
554289bc588SBarry Smith   else if (type == NORM_1) {
555289bc588SBarry Smith     *norm = 0.0;
556289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
557289bc588SBarry Smith       sum = 0.0;
558289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
559289bc588SBarry Smith #if defined(PETSC_COMPLEX)
560289bc588SBarry Smith         sum += abs(*v++);
561289bc588SBarry Smith #else
562289bc588SBarry Smith         sum += fabs(*v++);
563289bc588SBarry Smith #endif
564289bc588SBarry Smith       }
565289bc588SBarry Smith       if (sum > *norm) *norm = sum;
566289bc588SBarry Smith     }
567289bc588SBarry Smith   }
568289bc588SBarry Smith   else if (type == NORM_INFINITY) {
569289bc588SBarry Smith     *norm = 0.0;
570289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
571289bc588SBarry Smith       v = mat->v + j;
572289bc588SBarry Smith       sum = 0.0;
573289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
574289bc588SBarry Smith #if defined(PETSC_COMPLEX)
575289bc588SBarry Smith         sum += abs(*v); v += mat->m;
576289bc588SBarry Smith #else
577289bc588SBarry Smith         sum += fabs(*v); v += mat->m;
578289bc588SBarry Smith #endif
579289bc588SBarry Smith       }
580289bc588SBarry Smith       if (sum > *norm) *norm = sum;
581289bc588SBarry Smith     }
582289bc588SBarry Smith   }
583289bc588SBarry Smith   else {
58448d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
585289bc588SBarry Smith   }
586289bc588SBarry Smith   return 0;
587289bc588SBarry Smith }
588289bc588SBarry Smith 
589c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
590289bc588SBarry Smith {
591c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
592289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
593289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
594c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
595c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
596c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
597c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
598c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
599c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
600c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
601c0bbcb79SLois Curfman McInnes     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n");
602c0bbcb79SLois Curfman McInnes   else
6034d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
604289bc588SBarry Smith   return 0;
605289bc588SBarry Smith }
606289bc588SBarry Smith 
607ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
6086f0a148fSBarry Smith {
609ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
610416022c9SBarry Smith   PetscZero(l->v,l->m*l->n*sizeof(Scalar));
6116f0a148fSBarry Smith   return 0;
6126f0a148fSBarry Smith }
6136f0a148fSBarry Smith 
614ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
6156f0a148fSBarry Smith {
616ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
6176abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
6186f0a148fSBarry Smith   Scalar       *slot;
61978b31e54SBarry Smith   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
62078b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
6216f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
6226f0a148fSBarry Smith     slot = l->v + rows[i];
6236f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
6246f0a148fSBarry Smith   }
6256f0a148fSBarry Smith   if (diag) {
6266f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
6276f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
628260925b8SBarry Smith       *slot = *diag;
6296f0a148fSBarry Smith     }
6306f0a148fSBarry Smith   }
631260925b8SBarry Smith   ISRestoreIndices(is,&rows);
6326f0a148fSBarry Smith   return 0;
6336f0a148fSBarry Smith }
634557bce09SLois Curfman McInnes 
635c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
636557bce09SLois Curfman McInnes {
637c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
638557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
639557bce09SLois Curfman McInnes   return 0;
640557bce09SLois Curfman McInnes }
641557bce09SLois Curfman McInnes 
642c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
643d3e5ee88SLois Curfman McInnes {
644c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
645d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
646d3e5ee88SLois Curfman McInnes   return 0;
647d3e5ee88SLois Curfman McInnes }
648d3e5ee88SLois Curfman McInnes 
649c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
65064e87e97SBarry Smith {
651c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
65264e87e97SBarry Smith   *array = mat->v;
65364e87e97SBarry Smith   return 0;
65464e87e97SBarry Smith }
6550754003eSLois Curfman McInnes 
6560754003eSLois Curfman McInnes 
657c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
6580754003eSLois Curfman McInnes {
659ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
6600754003eSLois Curfman McInnes }
6610754003eSLois Curfman McInnes 
662c0bbcb79SLois Curfman McInnes static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,Mat *submat)
6630754003eSLois Curfman McInnes {
664c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
6650754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
666160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
6670754003eSLois Curfman McInnes   Scalar       *vwork, *val;
6680754003eSLois Curfman McInnes   Mat          newmat;
6690754003eSLois Curfman McInnes 
67078b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
67178b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
67278b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
67378b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
6740754003eSLois Curfman McInnes 
67578b31e54SBarry Smith   smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap);
67678b31e54SBarry Smith   cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork);
67778b31e54SBarry Smith   vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
678416022c9SBarry Smith   PetscZero((char*)smap,oldcols*sizeof(int));
6790754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
6800754003eSLois Curfman McInnes 
6810754003eSLois Curfman McInnes   /* Create and fill new matrix */
682c0bbcb79SLois Curfman McInnes   ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr);
6830754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
6840754003eSLois Curfman McInnes     nznew = 0;
6850754003eSLois Curfman McInnes     val   = mat->v + irow[i];
6860754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
6870754003eSLois Curfman McInnes       if (smap[j]) {
6880754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
6890754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
6900754003eSLois Curfman McInnes       }
6910754003eSLois Curfman McInnes     }
692dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
69378b31e54SBarry Smith            CHKERRQ(ierr);
6940754003eSLois Curfman McInnes   }
69578b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
69678b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
6970754003eSLois Curfman McInnes 
6980754003eSLois Curfman McInnes   /* Free work space */
69978b31e54SBarry Smith   PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork);
70078b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
70178b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
7020754003eSLois Curfman McInnes   *submat = newmat;
7030754003eSLois Curfman McInnes   return 0;
7040754003eSLois Curfman McInnes }
7050754003eSLois Curfman McInnes 
706289bc588SBarry Smith /* -------------------------------------------------------------------*/
707ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense,
708ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
709ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
710ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
711ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
712ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
713ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
714ec8511deSBarry Smith        MatRelax_SeqDense,
715ec8511deSBarry Smith        MatTranspose_SeqDense,
716ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
717ec8511deSBarry Smith        MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense,
718289bc588SBarry Smith        0,0,
719ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
720ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
721ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
722ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
723ec8511deSBarry Smith        0,0,MatGetArray_SeqDense,0,0,
724ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
7251987afe7SBarry Smith        MatCopyPrivate_SeqDense,0,0,0,0,
7261987afe7SBarry Smith        MatAXPY_SeqDense};
7270754003eSLois Curfman McInnes 
7284b828684SBarry Smith /*@C
729fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
730d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
731d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
732289bc588SBarry Smith 
73320563c6bSBarry Smith    Input Parameters:
7340c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
7350c775827SLois Curfman McInnes .  m - number of rows
7360c775827SLois Curfman McInnes .  n - number of column
73720563c6bSBarry Smith 
73820563c6bSBarry Smith    Output Parameter:
7390c775827SLois Curfman McInnes .  newmat - the matrix
74020563c6bSBarry Smith 
741dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
742d65003e9SLois Curfman McInnes 
743dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
74420563c6bSBarry Smith @*/
745fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat)
746289bc588SBarry Smith {
747ec8511deSBarry Smith   int          size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar);
748289bc588SBarry Smith   Mat          mat;
749ec8511deSBarry Smith   Mat_SeqDense *l;
750289bc588SBarry Smith   *newmat        = 0;
751ec8511deSBarry Smith   PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
752a5a9c739SBarry Smith   PLogObjectCreate(mat);
753ec8511deSBarry Smith   l              = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l);
754416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
755ec8511deSBarry Smith   mat->destroy   = MatDestroy_SeqDense;
756ec8511deSBarry Smith   mat->view      = MatView_SeqDense;
757289bc588SBarry Smith   mat->data      = (void *) l;
758289bc588SBarry Smith   mat->factor    = 0;
759752f5567SLois Curfman McInnes   PLogObjectMemory(mat,sizeof(struct _Mat) + size);
760d6dfbf8fSBarry Smith 
761289bc588SBarry Smith   l->m           = m;
762289bc588SBarry Smith   l->n           = n;
763289bc588SBarry Smith   l->v           = (Scalar *) (l + 1);
764289bc588SBarry Smith   l->pivots      = 0;
765289bc588SBarry Smith   l->roworiented = 1;
766d6dfbf8fSBarry Smith 
767416022c9SBarry Smith   PetscZero(l->v,m*n*sizeof(Scalar));
768289bc588SBarry Smith   *newmat = mat;
769289bc588SBarry Smith   return 0;
770289bc588SBarry Smith }
771289bc588SBarry Smith 
772c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
773289bc588SBarry Smith {
774c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
775c0bbcb79SLois Curfman McInnes   return MatCreateSeqDense(A->comm,m->m,m->n,newmat);
776289bc588SBarry Smith }
777