xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 94a424c13da5ee4d79d354b02fcfacbd9f84bab2)
1cb512458SBarry Smith #ifndef lint
2*94a424c1SBarry Smith static char vcid[] = "$Id: dense.c,v 1.98 1996/03/21 00:44:51 bsmith Exp bsmith $";
3cb512458SBarry Smith #endif
467e560aaSBarry Smith /*
567e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
667e560aaSBarry Smith */
7289bc588SBarry Smith 
817699dbbSLois Curfman McInnes #include "dense.h"
9b16a3bb1SBarry Smith #include "pinclude/plapack.h"
10b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
11289bc588SBarry Smith 
121987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data;
151987afe7SBarry Smith   int          N = x->m*x->n, one = 1;
161987afe7SBarry Smith   BLaxpy_( &N, alpha, x->v, &one, y->v, &one );
171987afe7SBarry Smith   PLogFlops(2*N-1);
181987afe7SBarry Smith   return 0;
191987afe7SBarry Smith }
201987afe7SBarry Smith 
21c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
22289bc588SBarry Smith {
23c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
24289bc588SBarry Smith   int          i,N = mat->m*mat->n,count = 0;
25289bc588SBarry Smith   Scalar       *v = mat->v;
26289bc588SBarry Smith   for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;}
27bcd2baecSBarry Smith   if (nz)      *nz      = count;
28bcd2baecSBarry Smith   if (nzalloc) *nzalloc = N;
29bcd2baecSBarry Smith   if (mem)     *mem     = (int)A->mem;
30fa9ec3f1SLois Curfman McInnes   return 0;
31289bc588SBarry Smith }
32289bc588SBarry Smith 
33289bc588SBarry Smith /* ---------------------------------------------------------------*/
34289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots,
35289bc588SBarry Smith    rather than put it in the Mat->col slot.*/
36c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f)
37289bc588SBarry Smith {
38c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
396abc6512SBarry Smith   int          info;
4055659b69SBarry Smith 
41289bc588SBarry Smith   if (!mat->pivots) {
420452661fSBarry Smith     mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots);
43c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,mat->m*sizeof(int));
44289bc588SBarry Smith   }
45289bc588SBarry Smith   LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info);
46ec8511deSBarry Smith   if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization");
47c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_LU;
4855659b69SBarry Smith   PLogFlops((2*mat->n*mat->n*mat->n)/3);
49289bc588SBarry Smith   return 0;
50289bc588SBarry Smith }
51c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact)
52289bc588SBarry Smith {
53a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
54289bc588SBarry Smith }
55c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
56289bc588SBarry Smith {
5749d8b64dSBarry Smith   return MatLUFactor(*fact,0,0,1.0);
58289bc588SBarry Smith }
59c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact)
60289bc588SBarry Smith {
61a57ad546SLois Curfman McInnes   return MatConvert(A,MATSAME,fact);
62289bc588SBarry Smith }
63c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
64289bc588SBarry Smith {
6549d8b64dSBarry Smith   return MatCholeskyFactor(*fact,0,1.0);
66289bc588SBarry Smith }
67c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f)
68289bc588SBarry Smith {
69c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
706abc6512SBarry Smith   int           info;
7155659b69SBarry Smith 
72752f5567SLois Curfman McInnes   if (mat->pivots) {
730452661fSBarry Smith     PetscFree(mat->pivots);
74c0bbcb79SLois Curfman McInnes     PLogObjectMemory(A,-mat->m*sizeof(int));
75752f5567SLois Curfman McInnes     mat->pivots = 0;
76752f5567SLois Curfman McInnes   }
77289bc588SBarry Smith   LApotrf_("L",&mat->n,mat->v,&mat->m,&info);
78ec8511deSBarry Smith   if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization");
79c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
8055659b69SBarry Smith   PLogFlops((mat->n*mat->n*mat->n)/3);
81289bc588SBarry Smith   return 0;
82289bc588SBarry Smith }
83289bc588SBarry Smith 
84c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
85289bc588SBarry Smith {
86c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
87a57ad546SLois Curfman McInnes   int          one = 1, info, ierr;
886abc6512SBarry Smith   Scalar       *x, *y;
8967e560aaSBarry Smith 
90a57ad546SLois Curfman McInnes   ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr);
91416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
92c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
9348d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
94289bc588SBarry Smith   }
95c0bbcb79SLois Curfman McInnes   else if (A->factor == FACTOR_CHOLESKY){
9648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
97289bc588SBarry Smith   }
98ec8511deSBarry Smith   else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve");
99ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve");
10055659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
101289bc588SBarry Smith   return 0;
102289bc588SBarry Smith }
103c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy)
104da3a660dSBarry Smith {
105c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1066abc6512SBarry Smith   int          one = 1, info;
1076abc6512SBarry Smith   Scalar       *x, *y;
10867e560aaSBarry Smith 
109da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
110416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
111752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
112da3a660dSBarry Smith   if (mat->pivots) {
11348d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
114da3a660dSBarry Smith   }
115da3a660dSBarry Smith   else {
11648d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
117da3a660dSBarry Smith   }
118ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve");
11955659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
120da3a660dSBarry Smith   return 0;
121da3a660dSBarry Smith }
122c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
123da3a660dSBarry Smith {
124c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
1256abc6512SBarry Smith   int          one = 1, info,ierr;
1266abc6512SBarry Smith   Scalar       *x, *y, sone = 1.0;
127da3a660dSBarry Smith   Vec          tmp = 0;
12867e560aaSBarry Smith 
129da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
130da3a660dSBarry Smith   if (yy == zz) {
13178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
132c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
13378b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
134da3a660dSBarry Smith   }
135416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
136752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
137da3a660dSBarry Smith   if (mat->pivots) {
13848d91487SBarry Smith     LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
139da3a660dSBarry Smith   }
140da3a660dSBarry Smith   else {
14148d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
142da3a660dSBarry Smith   }
143ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve");
144da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
145da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
14655659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
147da3a660dSBarry Smith   return 0;
148da3a660dSBarry Smith }
14967e560aaSBarry Smith 
150c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy)
151da3a660dSBarry Smith {
152c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense *) A->data;
1536abc6512SBarry Smith   int           one = 1, info,ierr;
1546abc6512SBarry Smith   Scalar        *x, *y, sone = 1.0;
155da3a660dSBarry Smith   Vec           tmp;
15667e560aaSBarry Smith 
157da3a660dSBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
158da3a660dSBarry Smith   if (yy == zz) {
15978b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr);
160c0bbcb79SLois Curfman McInnes     PLogObjectParent(A,tmp);
16178b31e54SBarry Smith     ierr = VecCopy(yy,tmp); CHKERRQ(ierr);
162da3a660dSBarry Smith   }
163416022c9SBarry Smith   PetscMemcpy(y,x,mat->m*sizeof(Scalar));
164752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
165da3a660dSBarry Smith   if (mat->pivots) {
16648d91487SBarry Smith     LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info );
167da3a660dSBarry Smith   }
168da3a660dSBarry Smith   else {
16948d91487SBarry Smith     LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info );
170da3a660dSBarry Smith   }
171ec8511deSBarry Smith   if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve");
172da3a660dSBarry Smith   if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);}
173da3a660dSBarry Smith   else VecAXPY(&sone,zz,yy);
17455659b69SBarry Smith   PLogFlops(mat->n*mat->n - mat->n);
175da3a660dSBarry Smith   return 0;
176da3a660dSBarry Smith }
177289bc588SBarry Smith /* ------------------------------------------------------------------*/
178c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag,
17920e1a0d4SLois Curfman McInnes                           double shift,int its,Vec xx)
180289bc588SBarry Smith {
181c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
182289bc588SBarry Smith   Scalar       *x, *b, *v = mat->v, zero = 0.0, xt;
1836abc6512SBarry Smith   int          o = 1,ierr, m = mat->m, i;
184289bc588SBarry Smith 
185289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
186289bc588SBarry Smith     /* this is a hack fix, should have another version without
187289bc588SBarry Smith        the second BLdot */
188bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx); CHKERRQ(ierr);
189289bc588SBarry Smith   }
190289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(bb,&b);
191289bc588SBarry Smith   while (its--) {
192289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
193289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
194f1747703SBarry Smith #if defined(PETSC_COMPLEX)
195f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
196f1747703SBarry Smith            not happy about returning a double complex */
197f1747703SBarry Smith         int    _i;
198f1747703SBarry Smith         Scalar sum = b[i];
199f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
200f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
201f1747703SBarry Smith         }
202f1747703SBarry Smith         xt = sum;
203f1747703SBarry Smith #else
204289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
205f1747703SBarry Smith #endif
206d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
207289bc588SBarry Smith       }
208289bc588SBarry Smith     }
209289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
210289bc588SBarry Smith       for ( i=m-1; i>=0; i-- ) {
211f1747703SBarry Smith #if defined(PETSC_COMPLEX)
212f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
213f1747703SBarry Smith            not happy about returning a double complex */
214f1747703SBarry Smith         int    _i;
215f1747703SBarry Smith         Scalar sum = b[i];
216f1747703SBarry Smith         for ( _i=0; _i<m; _i++ ) {
217f1747703SBarry Smith           sum -= conj(v[i+_i*m])*x[_i];
218f1747703SBarry Smith         }
219f1747703SBarry Smith         xt = sum;
220f1747703SBarry Smith #else
221289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
222f1747703SBarry Smith #endif
223d6dfbf8fSBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]);
224289bc588SBarry Smith       }
225289bc588SBarry Smith     }
226289bc588SBarry Smith   }
227289bc588SBarry Smith   return 0;
228289bc588SBarry Smith }
229289bc588SBarry Smith 
230289bc588SBarry Smith /* -----------------------------------------------------------------*/
231c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy)
232289bc588SBarry Smith {
233c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
234289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
235289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
236289bc588SBarry Smith   VecGetArray(xx,&x), VecGetArray(yy,&y);
23748d91487SBarry Smith   LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One);
23855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->n);
239289bc588SBarry Smith   return 0;
240289bc588SBarry Smith }
241c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
242289bc588SBarry Smith {
243c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
244289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y;
245289bc588SBarry Smith   int          _One=1;Scalar _DOne=1.0, _DZero=0.0;
246289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
24748d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One);
24855659b69SBarry Smith   PLogFlops(2*mat->m*mat->n - mat->m);
249289bc588SBarry Smith   return 0;
250289bc588SBarry Smith }
251c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
252289bc588SBarry Smith {
253c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
254289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
2556abc6512SBarry Smith   int          _One=1; Scalar _DOne=1.0;
256289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
257416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
25848d91487SBarry Smith   LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One);
25955659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
260289bc588SBarry Smith   return 0;
261289bc588SBarry Smith }
262c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
263289bc588SBarry Smith {
264c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
265289bc588SBarry Smith   Scalar       *v = mat->v, *x, *y, *z;
266289bc588SBarry Smith   int          _One=1;
2676abc6512SBarry Smith   Scalar       _DOne=1.0;
268289bc588SBarry Smith   VecGetArray(xx,&x); VecGetArray(yy,&y);
269289bc588SBarry Smith   VecGetArray(zz,&z);
270416022c9SBarry Smith   if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar));
27148d91487SBarry Smith   LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One);
27255659b69SBarry Smith   PLogFlops(2*mat->m*mat->n);
273289bc588SBarry Smith   return 0;
274289bc588SBarry Smith }
275289bc588SBarry Smith 
276289bc588SBarry Smith /* -----------------------------------------------------------------*/
277c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
278289bc588SBarry Smith {
279c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
280289bc588SBarry Smith   Scalar       *v;
281289bc588SBarry Smith   int          i;
28267e560aaSBarry Smith 
283289bc588SBarry Smith   *ncols = mat->n;
284289bc588SBarry Smith   if (cols) {
2850452661fSBarry Smith     *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols);
28673c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) (*cols)[i] = i;
287289bc588SBarry Smith   }
288289bc588SBarry Smith   if (vals) {
2890452661fSBarry Smith     *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals);
290289bc588SBarry Smith     v = mat->v + row;
29173c3ce57SBarry Smith     for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;}
292289bc588SBarry Smith   }
293289bc588SBarry Smith   return 0;
294289bc588SBarry Smith }
295c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals)
296289bc588SBarry Smith {
2970452661fSBarry Smith   if (cols) { PetscFree(*cols); }
2980452661fSBarry Smith   if (vals) { PetscFree(*vals); }
299289bc588SBarry Smith   return 0;
300289bc588SBarry Smith }
301289bc588SBarry Smith /* ----------------------------------------------------------------*/
302ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n,
303e8d4e0b9SBarry Smith                                     int *indexn,Scalar *v,InsertMode addv)
304289bc588SBarry Smith {
305c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
306289bc588SBarry Smith   int          i,j;
307d6dfbf8fSBarry Smith 
308289bc588SBarry Smith   if (!mat->roworiented) {
309dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
310289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
311289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
312289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
313289bc588SBarry Smith         }
314289bc588SBarry Smith       }
315289bc588SBarry Smith     }
316289bc588SBarry Smith     else {
317289bc588SBarry Smith       for ( j=0; j<n; j++ ) {
318289bc588SBarry Smith         for ( i=0; i<m; i++ ) {
319289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
320289bc588SBarry Smith         }
321289bc588SBarry Smith       }
322289bc588SBarry Smith     }
323e8d4e0b9SBarry Smith   }
324e8d4e0b9SBarry Smith   else {
325dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
326e8d4e0b9SBarry Smith       for ( i=0; i<m; i++ ) {
327e8d4e0b9SBarry Smith         for ( j=0; j<n; j++ ) {
328e8d4e0b9SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] = *v++;
329e8d4e0b9SBarry Smith         }
330e8d4e0b9SBarry Smith       }
331e8d4e0b9SBarry Smith     }
332289bc588SBarry Smith     else {
333289bc588SBarry Smith       for ( i=0; i<m; i++ ) {
334289bc588SBarry Smith         for ( j=0; j<n; j++ ) {
335289bc588SBarry Smith           mat->v[indexn[j]*mat->m + indexm[i]] += *v++;
336289bc588SBarry Smith         }
337289bc588SBarry Smith       }
338289bc588SBarry Smith     }
339e8d4e0b9SBarry Smith   }
340289bc588SBarry Smith   return 0;
341289bc588SBarry Smith }
342e8d4e0b9SBarry Smith 
343ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v)
344ae80bb75SLois Curfman McInnes {
345ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
346ae80bb75SLois Curfman McInnes   int          i, j;
347ae80bb75SLois Curfman McInnes   Scalar       *vpt = v;
348ae80bb75SLois Curfman McInnes 
349ae80bb75SLois Curfman McInnes   /* row-oriented output */
350ae80bb75SLois Curfman McInnes   for ( i=0; i<m; i++ ) {
351ae80bb75SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
352ae80bb75SLois Curfman McInnes       *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]];
353ae80bb75SLois Curfman McInnes     }
354ae80bb75SLois Curfman McInnes   }
355ae80bb75SLois Curfman McInnes   return 0;
356ae80bb75SLois Curfman McInnes }
357ae80bb75SLois Curfman McInnes 
358289bc588SBarry Smith /* -----------------------------------------------------------------*/
3593d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues)
360289bc588SBarry Smith {
361c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l;
362289bc588SBarry Smith   int          ierr;
363289bc588SBarry Smith   Mat          newi;
36448d91487SBarry Smith 
365b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr);
366ec8511deSBarry Smith   l = (Mat_SeqDense *) newi->data;
36755659b69SBarry Smith   if (cpvalues == COPY_VALUES) {
368416022c9SBarry Smith     PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));
36955659b69SBarry Smith   }
370227d817aSBarry Smith   newi->assembled = PETSC_TRUE;
371289bc588SBarry Smith   *newmat = newi;
372289bc588SBarry Smith   return 0;
373289bc588SBarry Smith }
374289bc588SBarry Smith 
37577c4ece6SBarry Smith #include "sys.h"
37688e32edaSLois Curfman McInnes 
37719bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A)
37888e32edaSLois Curfman McInnes {
37988e32edaSLois Curfman McInnes   Mat_SeqDense *a;
38088e32edaSLois Curfman McInnes   Mat          B;
38188e32edaSLois Curfman McInnes   int          *scols, i, j, nz, ierr, fd, header[4], size;
38288e32edaSLois Curfman McInnes   int          *rowlengths = 0, M, N, *cols;
38388e32edaSLois Curfman McInnes   Scalar       *vals, *svals, *v;
38419bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
38588e32edaSLois Curfman McInnes 
38688e32edaSLois Curfman McInnes   MPI_Comm_size(comm,&size);
38788e32edaSLois Curfman McInnes   if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor");
38890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
38977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
39088e32edaSLois Curfman McInnes   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object");
39188e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
39288e32edaSLois Curfman McInnes 
39390ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
39490ace30eSBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
39590ace30eSBarry Smith     B = *A;
39690ace30eSBarry Smith     a = (Mat_SeqDense *) B->data;
39790ace30eSBarry Smith 
39890ace30eSBarry Smith     /* read in nonzero values */
39977c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr);
40090ace30eSBarry Smith 
40190ace30eSBarry Smith     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
40290ace30eSBarry Smith     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
40390ace30eSBarry Smith     ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr);
40490ace30eSBarry Smith   } else {
40588e32edaSLois Curfman McInnes     /* read row lengths */
4060452661fSBarry Smith     rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
40777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
40888e32edaSLois Curfman McInnes 
40988e32edaSLois Curfman McInnes     /* create our matrix */
410b4fd4287SBarry Smith     ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr);
41188e32edaSLois Curfman McInnes     B = *A;
41288e32edaSLois Curfman McInnes     a = (Mat_SeqDense *) B->data;
41388e32edaSLois Curfman McInnes     v = a->v;
41488e32edaSLois Curfman McInnes 
41588e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
4160452661fSBarry Smith     cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols);
41777c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr);
4180452661fSBarry Smith     vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals);
41977c4ece6SBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr);
42088e32edaSLois Curfman McInnes 
42188e32edaSLois Curfman McInnes     /* insert into matrix */
42288e32edaSLois Curfman McInnes     for ( i=0; i<M; i++ ) {
42388e32edaSLois Curfman McInnes       for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j];
42488e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
42588e32edaSLois Curfman McInnes     }
4260452661fSBarry Smith     PetscFree(vals); PetscFree(cols); PetscFree(rowlengths);
42788e32edaSLois Curfman McInnes 
42888e32edaSLois Curfman McInnes     ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
42988e32edaSLois Curfman McInnes     ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
43090ace30eSBarry Smith   }
43188e32edaSLois Curfman McInnes   return 0;
43288e32edaSLois Curfman McInnes }
43388e32edaSLois Curfman McInnes 
434932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h"
43577c4ece6SBarry Smith #include "sys.h"
436932b0c3eSLois Curfman McInnes 
437932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer)
438289bc588SBarry Smith {
439932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
440932b0c3eSLois Curfman McInnes   int          ierr, i, j, format;
441d636dbe3SBarry Smith   FILE         *fd;
442932b0c3eSLois Curfman McInnes   char         *outputname;
443932b0c3eSLois Curfman McInnes   Scalar       *v;
444932b0c3eSLois Curfman McInnes 
44590ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
446932b0c3eSLois Curfman McInnes   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
44790ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
44890ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
449932b0c3eSLois Curfman McInnes     ;  /* do nothing for now */
450f72585f2SLois Curfman McInnes   }
451f72585f2SLois Curfman McInnes   else {
45247989497SBarry Smith #if defined(PETSC_COMPLEX)
45347989497SBarry Smith     int allreal = 1;
45447989497SBarry Smith     /* determine if matrix has all real values */
45547989497SBarry Smith     v = a->v;
45647989497SBarry Smith     for ( i=0; i<a->m*a->n; i++ ) {
45747989497SBarry Smith       if (imag(v[i])) { allreal = 0; break ;}
45847989497SBarry Smith     }
45947989497SBarry Smith #endif
460932b0c3eSLois Curfman McInnes     for ( i=0; i<a->m; i++ ) {
461932b0c3eSLois Curfman McInnes       v = a->v + i;
462932b0c3eSLois Curfman McInnes       for ( j=0; j<a->n; j++ ) {
463289bc588SBarry Smith #if defined(PETSC_COMPLEX)
46447989497SBarry Smith         if (allreal) {
46547989497SBarry Smith           fprintf(fd,"%6.4e ",real(*v)); v += a->m;
46647989497SBarry Smith         } else {
467932b0c3eSLois Curfman McInnes           fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m;
46847989497SBarry Smith         }
469289bc588SBarry Smith #else
470932b0c3eSLois Curfman McInnes         fprintf(fd,"%6.4e ",*v); v += a->m;
471289bc588SBarry Smith #endif
472289bc588SBarry Smith       }
4738e37dc09SBarry Smith       fprintf(fd,"\n");
474289bc588SBarry Smith     }
475da3a660dSBarry Smith   }
476932b0c3eSLois Curfman McInnes   fflush(fd);
477289bc588SBarry Smith   return 0;
478289bc588SBarry Smith }
479289bc588SBarry Smith 
480932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer)
481932b0c3eSLois Curfman McInnes {
482932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense *) A->data;
483932b0c3eSLois Curfman McInnes   int          ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n;
48490ace30eSBarry Smith   int          format;
48590ace30eSBarry Smith   Scalar       *v, *anonz,*vals;
486932b0c3eSLois Curfman McInnes 
48790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
48890ace30eSBarry Smith 
48990ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
49090ace30eSBarry Smith   if (format == BINARY_FORMAT_NATIVE) {
49190ace30eSBarry Smith     /* store the matrix as a dense matrix */
49290ace30eSBarry Smith     col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens);
49390ace30eSBarry Smith     col_lens[0] = MAT_COOKIE;
49490ace30eSBarry Smith     col_lens[1] = m;
49590ace30eSBarry Smith     col_lens[2] = n;
49690ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
49777c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr);
49890ace30eSBarry Smith     PetscFree(col_lens);
49990ace30eSBarry Smith 
50090ace30eSBarry Smith     /* write out matrix, by rows */
50190ace30eSBarry Smith     vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals);
50290ace30eSBarry Smith     v    = a->v;
50390ace30eSBarry Smith     for ( i=0; i<m; i++ ) {
50490ace30eSBarry Smith       for ( j=0; j<n; j++ ) {
50590ace30eSBarry Smith         vals[i + j*m] = *v++;
50690ace30eSBarry Smith       }
50790ace30eSBarry Smith     }
50877c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr);
50990ace30eSBarry Smith     PetscFree(vals);
51090ace30eSBarry Smith   } else {
5110452661fSBarry Smith     col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens);
512932b0c3eSLois Curfman McInnes     col_lens[0] = MAT_COOKIE;
513932b0c3eSLois Curfman McInnes     col_lens[1] = m;
514932b0c3eSLois Curfman McInnes     col_lens[2] = n;
515932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
516932b0c3eSLois Curfman McInnes 
517932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
518932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) col_lens[4+i] = n;
51977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr);
520932b0c3eSLois Curfman McInnes 
521932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
522932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
523932b0c3eSLois Curfman McInnes     ict = 0;
524932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
525932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) col_lens[ict++] = j;
526932b0c3eSLois Curfman McInnes     }
52777c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr);
5280452661fSBarry Smith     PetscFree(col_lens);
529932b0c3eSLois Curfman McInnes 
530932b0c3eSLois Curfman McInnes     /* store nonzero values */
5310452661fSBarry Smith     anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz);
532932b0c3eSLois Curfman McInnes     ict = 0;
533932b0c3eSLois Curfman McInnes     for ( i=0; i<m; i++ ) {
534932b0c3eSLois Curfman McInnes       v = a->v + i;
535932b0c3eSLois Curfman McInnes       for ( j=0; j<n; j++ ) {
536932b0c3eSLois Curfman McInnes         anonz[ict++] = *v; v += a->m;
537932b0c3eSLois Curfman McInnes       }
538932b0c3eSLois Curfman McInnes     }
53977c4ece6SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr);
5400452661fSBarry Smith     PetscFree(anonz);
54190ace30eSBarry Smith   }
542932b0c3eSLois Curfman McInnes   return 0;
543932b0c3eSLois Curfman McInnes }
544932b0c3eSLois Curfman McInnes 
545932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer)
546932b0c3eSLois Curfman McInnes {
547932b0c3eSLois Curfman McInnes   Mat          A = (Mat) obj;
548932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*) A->data;
549bcd2baecSBarry Smith   ViewerType   vtype;
550bcd2baecSBarry Smith   int          ierr;
551932b0c3eSLois Curfman McInnes 
552932b0c3eSLois Curfman McInnes   if (!viewer) {
55319bcc07fSBarry Smith     viewer = STDOUT_VIEWER_SELF;
554932b0c3eSLois Curfman McInnes   }
555bcd2baecSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
556bcd2baecSBarry Smith 
557bcd2baecSBarry Smith   if (vtype == MATLAB_VIEWER) {
558932b0c3eSLois Curfman McInnes     return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v);
559932b0c3eSLois Curfman McInnes   }
56019bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
561932b0c3eSLois Curfman McInnes     return MatView_SeqDense_ASCII(A,viewer);
562932b0c3eSLois Curfman McInnes   }
563bcd2baecSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
564932b0c3eSLois Curfman McInnes     return MatView_SeqDense_Binary(A,viewer);
565932b0c3eSLois Curfman McInnes   }
566932b0c3eSLois Curfman McInnes   return 0;
567932b0c3eSLois Curfman McInnes }
568289bc588SBarry Smith 
569ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj)
570289bc588SBarry Smith {
571289bc588SBarry Smith   Mat          mat = (Mat) obj;
572ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) mat->data;
573a5a9c739SBarry Smith #if defined(PETSC_LOG)
574a5a9c739SBarry Smith   PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n);
575a5a9c739SBarry Smith #endif
5760452661fSBarry Smith   if (l->pivots) PetscFree(l->pivots);
5773439631bSBarry Smith   if (!l->user_alloc) PetscFree(l->v);
5780452661fSBarry Smith   PetscFree(l);
579a5a9c739SBarry Smith   PLogObjectDestroy(mat);
5800452661fSBarry Smith   PetscHeaderDestroy(mat);
581289bc588SBarry Smith   return 0;
582289bc588SBarry Smith }
583289bc588SBarry Smith 
584c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout)
585289bc588SBarry Smith {
586c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
587d3e5ee88SLois Curfman McInnes   int          k, j, m, n;
588d3e5ee88SLois Curfman McInnes   Scalar       *v, tmp;
58948b35521SBarry Smith 
590d3e5ee88SLois Curfman McInnes   v = mat->v; m = mat->m; n = mat->n;
5913638b69dSLois Curfman McInnes   if (matout == PETSC_NULL) { /* in place transpose */
592d9f96c7cSLois Curfman McInnes     if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place");
593d3e5ee88SLois Curfman McInnes     for ( j=0; j<m; j++ ) {
594289bc588SBarry Smith       for ( k=0; k<j; k++ ) {
595d3e5ee88SLois Curfman McInnes         tmp = v[j + k*n];
596d3e5ee88SLois Curfman McInnes         v[j + k*n] = v[k + j*n];
597d3e5ee88SLois Curfman McInnes         v[k + j*n] = tmp;
598289bc588SBarry Smith       }
599289bc588SBarry Smith     }
60048b35521SBarry Smith   }
601d3e5ee88SLois Curfman McInnes   else { /* out-of-place transpose */
602d3e5ee88SLois Curfman McInnes     int          ierr;
603d3e5ee88SLois Curfman McInnes     Mat          tmat;
604ec8511deSBarry Smith     Mat_SeqDense *tmatd;
605d3e5ee88SLois Curfman McInnes     Scalar       *v2;
606b4fd4287SBarry Smith     ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr);
607ec8511deSBarry Smith     tmatd = (Mat_SeqDense *) tmat->data;
6080de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
609d3e5ee88SLois Curfman McInnes     for ( j=0; j<n; j++ ) {
6100de55854SLois Curfman McInnes       for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m];
611d3e5ee88SLois Curfman McInnes     }
612d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
613d3e5ee88SLois Curfman McInnes     ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
614d3e5ee88SLois Curfman McInnes     *matout = tmat;
61548b35521SBarry Smith   }
616289bc588SBarry Smith   return 0;
617289bc588SBarry Smith }
618289bc588SBarry Smith 
61977c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg)
620289bc588SBarry Smith {
621c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data;
622c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data;
623289bc588SBarry Smith   int          i;
624289bc588SBarry Smith   Scalar       *v1 = mat1->v, *v2 = mat2->v;
6259ea5d5aeSSatish Balay 
6269ea5d5aeSSatish Balay   if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type  MATSEQDENSE");
62777c4ece6SBarry Smith   if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;}
62877c4ece6SBarry Smith   if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;}
629289bc588SBarry Smith   for ( i=0; i<mat1->m*mat1->n; i++ ) {
63077c4ece6SBarry Smith     if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;}
631289bc588SBarry Smith     v1++; v2++;
632289bc588SBarry Smith   }
63377c4ece6SBarry Smith   *flg = PETSC_TRUE;
6349ea5d5aeSSatish Balay   return 0;
635289bc588SBarry Smith }
636289bc588SBarry Smith 
637c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v)
638289bc588SBarry Smith {
639c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
6406abc6512SBarry Smith   int          i, n;
6416abc6512SBarry Smith   Scalar       *x;
642289bc588SBarry Smith   VecGetArray(v,&x); VecGetSize(v,&n);
643ec8511deSBarry Smith   if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec");
644289bc588SBarry Smith   for ( i=0; i<mat->m; i++ ) {
645289bc588SBarry Smith     x[i] = mat->v[i*mat->m + i];
646289bc588SBarry Smith   }
647289bc588SBarry Smith   return 0;
648289bc588SBarry Smith }
649289bc588SBarry Smith 
650052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
651289bc588SBarry Smith {
652c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
653da3a660dSBarry Smith   Scalar       *l,*r,x,*v;
654da3a660dSBarry Smith   int          i,j,m = mat->m, n = mat->n;
65555659b69SBarry Smith 
65628988994SBarry Smith   if (ll) {
657da3a660dSBarry Smith     VecGetArray(ll,&l); VecGetSize(ll,&m);
658052efed2SBarry Smith     if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size");
65955659b69SBarry Smith     PLogFlops(n*m);
660da3a660dSBarry Smith     for ( i=0; i<m; i++ ) {
661da3a660dSBarry Smith       x = l[i];
662da3a660dSBarry Smith       v = mat->v + i;
663da3a660dSBarry Smith       for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;}
664da3a660dSBarry Smith     }
665da3a660dSBarry Smith   }
66628988994SBarry Smith   if (rr) {
667da3a660dSBarry Smith     VecGetArray(rr,&r); VecGetSize(rr,&n);
668052efed2SBarry Smith     if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size");
66955659b69SBarry Smith     PLogFlops(n*m);
670da3a660dSBarry Smith     for ( i=0; i<n; i++ ) {
671da3a660dSBarry Smith       x = r[i];
672da3a660dSBarry Smith       v = mat->v + i*m;
673da3a660dSBarry Smith       for ( j=0; j<m; j++ ) { (*v++) *= x;}
674da3a660dSBarry Smith     }
675da3a660dSBarry Smith   }
676289bc588SBarry Smith   return 0;
677289bc588SBarry Smith }
678289bc588SBarry Smith 
679cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm)
680289bc588SBarry Smith {
681c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
682289bc588SBarry Smith   Scalar       *v = mat->v;
683289bc588SBarry Smith   double       sum = 0.0;
684289bc588SBarry Smith   int          i, j;
68555659b69SBarry Smith 
686289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
687289bc588SBarry Smith     for (i=0; i<mat->n*mat->m; i++ ) {
688289bc588SBarry Smith #if defined(PETSC_COMPLEX)
689289bc588SBarry Smith       sum += real(conj(*v)*(*v)); v++;
690289bc588SBarry Smith #else
691289bc588SBarry Smith       sum += (*v)*(*v); v++;
692289bc588SBarry Smith #endif
693289bc588SBarry Smith     }
694289bc588SBarry Smith     *norm = sqrt(sum);
69555659b69SBarry Smith     PLogFlops(2*mat->n*mat->m);
696289bc588SBarry Smith   }
697289bc588SBarry Smith   else if (type == NORM_1) {
698289bc588SBarry Smith     *norm = 0.0;
699289bc588SBarry Smith     for ( j=0; j<mat->n; j++ ) {
700289bc588SBarry Smith       sum = 0.0;
701289bc588SBarry Smith       for ( i=0; i<mat->m; i++ ) {
70233a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
703289bc588SBarry Smith       }
704289bc588SBarry Smith       if (sum > *norm) *norm = sum;
705289bc588SBarry Smith     }
70655659b69SBarry Smith     PLogFlops(mat->n*mat->m);
707289bc588SBarry Smith   }
708289bc588SBarry Smith   else if (type == NORM_INFINITY) {
709289bc588SBarry Smith     *norm = 0.0;
710289bc588SBarry Smith     for ( j=0; j<mat->m; j++ ) {
711289bc588SBarry Smith       v = mat->v + j;
712289bc588SBarry Smith       sum = 0.0;
713289bc588SBarry Smith       for ( i=0; i<mat->n; i++ ) {
714cddf8d76SBarry Smith         sum += PetscAbsScalar(*v); v += mat->m;
715289bc588SBarry Smith       }
716289bc588SBarry Smith       if (sum > *norm) *norm = sum;
717289bc588SBarry Smith     }
71855659b69SBarry Smith     PLogFlops(mat->n*mat->m);
719289bc588SBarry Smith   }
720289bc588SBarry Smith   else {
72148d91487SBarry Smith     SETERRQ(1,"MatNorm_SeqDense:No two norm");
722289bc588SBarry Smith   }
723289bc588SBarry Smith   return 0;
724289bc588SBarry Smith }
725289bc588SBarry Smith 
726c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op)
727289bc588SBarry Smith {
728c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense *) A->data;
72967e560aaSBarry Smith 
730289bc588SBarry Smith   if (op == ROW_ORIENTED)            aij->roworiented = 1;
731289bc588SBarry Smith   else if (op == COLUMN_ORIENTED)    aij->roworiented = 0;
732c0bbcb79SLois Curfman McInnes   else if (op == ROWS_SORTED ||
73388e32edaSLois Curfman McInnes            op == COLUMNS_SORTED ||
734c0bbcb79SLois Curfman McInnes            op == SYMMETRIC_MATRIX ||
735c0bbcb79SLois Curfman McInnes            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
736c0bbcb79SLois Curfman McInnes            op == NO_NEW_NONZERO_LOCATIONS ||
737c0bbcb79SLois Curfman McInnes            op == YES_NEW_NONZERO_LOCATIONS ||
738c0bbcb79SLois Curfman McInnes            op == NO_NEW_DIAGONALS ||
739c0bbcb79SLois Curfman McInnes            op == YES_NEW_DIAGONALS)
740*94a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n");
741c0bbcb79SLois Curfman McInnes   else
7424d39ef84SLois Curfman McInnes     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");}
743289bc588SBarry Smith   return 0;
744289bc588SBarry Smith }
745289bc588SBarry Smith 
746ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A)
7476f0a148fSBarry Smith {
748ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
749cddf8d76SBarry Smith   PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));
7506f0a148fSBarry Smith   return 0;
7516f0a148fSBarry Smith }
7526f0a148fSBarry Smith 
753ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag)
7546f0a148fSBarry Smith {
755ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense *) A->data;
7566abc6512SBarry Smith   int          n = l->n, i, j,ierr,N, *rows;
7576f0a148fSBarry Smith   Scalar       *slot;
75855659b69SBarry Smith 
75977c4ece6SBarry Smith   ierr = ISGetSize(is,&N); CHKERRQ(ierr);
76078b31e54SBarry Smith   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
7616f0a148fSBarry Smith   for ( i=0; i<N; i++ ) {
7626f0a148fSBarry Smith     slot = l->v + rows[i];
7636f0a148fSBarry Smith     for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;}
7646f0a148fSBarry Smith   }
7656f0a148fSBarry Smith   if (diag) {
7666f0a148fSBarry Smith     for ( i=0; i<N; i++ ) {
7676f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
768260925b8SBarry Smith       *slot = *diag;
7696f0a148fSBarry Smith     }
7706f0a148fSBarry Smith   }
771260925b8SBarry Smith   ISRestoreIndices(is,&rows);
7726f0a148fSBarry Smith   return 0;
7736f0a148fSBarry Smith }
774557bce09SLois Curfman McInnes 
775c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n)
776557bce09SLois Curfman McInnes {
777c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
778557bce09SLois Curfman McInnes   *m = mat->m; *n = mat->n;
779557bce09SLois Curfman McInnes   return 0;
780557bce09SLois Curfman McInnes }
781557bce09SLois Curfman McInnes 
782c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n)
783d3e5ee88SLois Curfman McInnes {
784c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
785d3e5ee88SLois Curfman McInnes   *m = 0; *n = mat->m;
786d3e5ee88SLois Curfman McInnes   return 0;
787d3e5ee88SLois Curfman McInnes }
788d3e5ee88SLois Curfman McInnes 
789c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array)
79064e87e97SBarry Smith {
791c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
79264e87e97SBarry Smith   *array = mat->v;
79364e87e97SBarry Smith   return 0;
79464e87e97SBarry Smith }
7950754003eSLois Curfman McInnes 
796ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array)
797ff14e315SSatish Balay {
798ff14e315SSatish Balay   return 0;
799ff14e315SSatish Balay }
8000754003eSLois Curfman McInnes 
801c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol)
8020754003eSLois Curfman McInnes {
803ec8511deSBarry Smith   SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done");
8040754003eSLois Curfman McInnes }
8050754003eSLois Curfman McInnes 
806cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,
807cddf8d76SBarry Smith                                     Mat *submat)
8080754003eSLois Curfman McInnes {
809c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense *) A->data;
8100754003eSLois Curfman McInnes   int          nznew, *smap, i, j, ierr, oldcols = mat->n;
811160018dcSBarry Smith   int          *irow, *icol, nrows, ncols, *cwork;
8120754003eSLois Curfman McInnes   Scalar       *vwork, *val;
8130754003eSLois Curfman McInnes   Mat          newmat;
8140754003eSLois Curfman McInnes 
81578b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
81678b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
81778b31e54SBarry Smith   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
81878b31e54SBarry Smith   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
8190754003eSLois Curfman McInnes 
8200452661fSBarry Smith   smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap);
8210452661fSBarry Smith   cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork);
8220452661fSBarry Smith   vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork);
823cddf8d76SBarry Smith   PetscMemzero((char*)smap,oldcols*sizeof(int));
8240754003eSLois Curfman McInnes   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
8250754003eSLois Curfman McInnes 
8260754003eSLois Curfman McInnes   /* Create and fill new matrix */
827b4fd4287SBarry Smith   ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr);
8280754003eSLois Curfman McInnes   for (i=0; i<nrows; i++) {
8290754003eSLois Curfman McInnes     nznew = 0;
8300754003eSLois Curfman McInnes     val   = mat->v + irow[i];
8310754003eSLois Curfman McInnes     for (j=0; j<oldcols; j++) {
8320754003eSLois Curfman McInnes       if (smap[j]) {
8330754003eSLois Curfman McInnes         cwork[nznew]   = smap[j] - 1;
8340754003eSLois Curfman McInnes         vwork[nznew++] = val[j * mat->m];
8350754003eSLois Curfman McInnes       }
8360754003eSLois Curfman McInnes     }
837dbb450caSBarry Smith     ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
83878b31e54SBarry Smith            CHKERRQ(ierr);
8390754003eSLois Curfman McInnes   }
84078b31e54SBarry Smith   ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
84178b31e54SBarry Smith   ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr);
8420754003eSLois Curfman McInnes 
8430754003eSLois Curfman McInnes   /* Free work space */
8440452661fSBarry Smith   PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
84578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
84678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
8470754003eSLois Curfman McInnes   *submat = newmat;
8480754003eSLois Curfman McInnes   return 0;
8490754003eSLois Curfman McInnes }
8500754003eSLois Curfman McInnes 
8514b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B)
8524b0e389bSBarry Smith {
8534b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data;
8544b0e389bSBarry Smith   if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B);
8554b0e389bSBarry Smith   if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)");
8564b0e389bSBarry Smith   PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));
8574b0e389bSBarry Smith   return 0;
8584b0e389bSBarry Smith }
8594b0e389bSBarry Smith 
860289bc588SBarry Smith /* -------------------------------------------------------------------*/
861ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense,
862ec8511deSBarry Smith        MatGetRow_SeqDense, MatRestoreRow_SeqDense,
863ec8511deSBarry Smith        MatMult_SeqDense, MatMultAdd_SeqDense,
864ec8511deSBarry Smith        MatMultTrans_SeqDense, MatMultTransAdd_SeqDense,
865ec8511deSBarry Smith        MatSolve_SeqDense,MatSolveAdd_SeqDense,
866ec8511deSBarry Smith        MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense,
867ec8511deSBarry Smith        MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense,
868ec8511deSBarry Smith        MatRelax_SeqDense,
869ec8511deSBarry Smith        MatTranspose_SeqDense,
870ec8511deSBarry Smith        MatGetInfo_SeqDense,MatEqual_SeqDense,
871052efed2SBarry Smith        MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense,
872289bc588SBarry Smith        0,0,
873ec8511deSBarry Smith        0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0,
874ec8511deSBarry Smith        MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense,
875ec8511deSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense,
876ec8511deSBarry Smith        MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense,
877ff14e315SSatish Balay        0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0,
878ec8511deSBarry Smith        MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense,
8793d1612f7SBarry Smith        MatConvertSameType_SeqDense,0,0,0,0,
880ae80bb75SLois Curfman McInnes        MatAXPY_SeqDense,0,0,
8814b0e389bSBarry Smith        MatGetValues_SeqDense,
8824b0e389bSBarry Smith        MatCopy_SeqDense};
8830754003eSLois Curfman McInnes 
88490ace30eSBarry Smith 
8854b828684SBarry Smith /*@C
886fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
887d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
888d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
889289bc588SBarry Smith 
89020563c6bSBarry Smith    Input Parameters:
8910c775827SLois Curfman McInnes .  comm - MPI communicator, set to MPI_COMM_SELF
8920c775827SLois Curfman McInnes .  m - number of rows
89318f449edSLois Curfman McInnes .  n - number of columns
894b4fd4287SBarry Smith .  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
895dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
89620563c6bSBarry Smith 
89720563c6bSBarry Smith    Output Parameter:
8980c775827SLois Curfman McInnes .  newmat - the matrix
89920563c6bSBarry Smith 
90018f449edSLois Curfman McInnes   Notes:
90118f449edSLois Curfman McInnes   The data input variable is intended primarily for Fortran programmers
90218f449edSLois Curfman McInnes   who wish to allocate their own matrix memory space.  Most users should
903b4fd4287SBarry Smith   set data=PETSC_NULL.
90418f449edSLois Curfman McInnes 
905dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
906d65003e9SLois Curfman McInnes 
907dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues()
90820563c6bSBarry Smith @*/
90918f449edSLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat)
910289bc588SBarry Smith {
911289bc588SBarry Smith   Mat          mat;
912ec8511deSBarry Smith   Mat_SeqDense *l;
91325cdf11fSBarry Smith   int          ierr,flg;
91455659b69SBarry Smith 
915289bc588SBarry Smith   *newmat        = 0;
91618f449edSLois Curfman McInnes 
9170452661fSBarry Smith   PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm);
918a5a9c739SBarry Smith   PLogObjectCreate(mat);
9193439631bSBarry Smith   l               = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l);
920416022c9SBarry Smith   PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps));
921ec8511deSBarry Smith   mat->destroy    = MatDestroy_SeqDense;
922ec8511deSBarry Smith   mat->view       = MatView_SeqDense;
923289bc588SBarry Smith   mat->factor     = 0;
9243439631bSBarry Smith   PLogObjectMemory(mat,sizeof(struct _Mat));
92518f449edSLois Curfman McInnes   mat->data       = (void *) l;
926d6dfbf8fSBarry Smith 
927289bc588SBarry Smith   l->m            = m;
928289bc588SBarry Smith   l->n            = n;
929289bc588SBarry Smith   l->pivots       = 0;
930289bc588SBarry Smith   l->roworiented  = 1;
931b4fd4287SBarry Smith   if (data == PETSC_NULL) {
932ba7af9b1SBarry Smith     l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v);
933cddf8d76SBarry Smith     PetscMemzero(l->v,m*n*sizeof(Scalar));
9342dd2b441SLois Curfman McInnes     l->user_alloc = 0;
9353439631bSBarry Smith     PLogObjectMemory(mat,n*m*sizeof(Scalar));
93618f449edSLois Curfman McInnes   }
9372dd2b441SLois Curfman McInnes   else { /* user-allocated storage */
9382dd2b441SLois Curfman McInnes     l->v = data;
9392dd2b441SLois Curfman McInnes     l->user_alloc = 1;
9402dd2b441SLois Curfman McInnes   }
94118f449edSLois Curfman McInnes 
94225cdf11fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
94325cdf11fSBarry Smith   if (flg) {
944682d7d0cSBarry Smith     ierr = MatPrintHelp(mat); CHKERRQ(ierr);
945682d7d0cSBarry Smith   }
946289bc588SBarry Smith   *newmat = mat;
947289bc588SBarry Smith   return 0;
948289bc588SBarry Smith }
949289bc588SBarry Smith 
950c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat)
951289bc588SBarry Smith {
952c0bbcb79SLois Curfman McInnes   Mat_SeqDense *m = (Mat_SeqDense *) A->data;
953b4fd4287SBarry Smith   return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat);
954289bc588SBarry Smith }
955