xref: /petsc/src/mat/impls/dense/seq/dense.c (revision f4df32b132e30afc79b77d8c3179eb18c28a5b40)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
367e560aaSBarry Smith /*
467e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
567e560aaSBarry Smith */
6289bc588SBarry Smith 
770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
8b4c862e0SSatish Balay #include "petscblaslapack.h"
9289bc588SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
12*f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
131987afe7SBarry Smith {
141987afe7SBarry Smith   Mat_SeqDense   *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
15*f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1613f74950SBarry Smith   PetscInt       j;
174ce68768SBarry Smith   PetscBLASInt   N = (PetscBLASInt)X->m*X->n,m=(PetscBLASInt)X->m,ldax = x->lda,lday=y->lda,one = 1;
18efee365bSSatish Balay   PetscErrorCode ierr;
193a40ed3dSBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
21*f4df32b1SMatthew Knepley   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(X) != size(Y)");
22a5ce6ee0Svictorle   if (ldax>m || lday>m) {
23a5ce6ee0Svictorle     for (j=0; j<X->n; j++) {
24*f4df32b1SMatthew Knepley       BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one);
25a5ce6ee0Svictorle     }
26a5ce6ee0Svictorle   } else {
27*f4df32b1SMatthew Knepley     BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one);
28a5ce6ee0Svictorle   }
29efee365bSSatish Balay   ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr);
303a40ed3dSBarry Smith   PetscFunctionReturn(0);
311987afe7SBarry Smith }
321987afe7SBarry Smith 
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
35dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
36289bc588SBarry Smith {
374e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3813f74950SBarry Smith   PetscInt     i,N = A->m*A->n,count = 0;
3987828ca2SBarry Smith   PetscScalar  *v = a->v;
403a40ed3dSBarry Smith 
413a40ed3dSBarry Smith   PetscFunctionBegin;
42289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
434e220ebcSLois Curfman McInnes 
44273d9f13SBarry Smith   info->rows_global       = (double)A->m;
45273d9f13SBarry Smith   info->columns_global    = (double)A->n;
46273d9f13SBarry Smith   info->rows_local        = (double)A->m;
47273d9f13SBarry Smith   info->columns_local     = (double)A->n;
484e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
494e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
504e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
514e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
524e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
534e220ebcSLois Curfman McInnes   info->mallocs           = 0;
544e220ebcSLois Curfman McInnes   info->memory            = A->mem;
554e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
564e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
574e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
584e220ebcSLois Curfman McInnes 
593a40ed3dSBarry Smith   PetscFunctionReturn(0);
60289bc588SBarry Smith }
61289bc588SBarry Smith 
624a2ae208SSatish Balay #undef __FUNCT__
634a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
64*f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
6580cd9d93SLois Curfman McInnes {
66273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
67*f4df32b1SMatthew Knepley   PetscScalar  oalpha = alpha;
684ce68768SBarry Smith   PetscBLASInt one = 1,lda = a->lda,j,nz;
69efee365bSSatish Balay   PetscErrorCode ierr;
7080cd9d93SLois Curfman McInnes 
713a40ed3dSBarry Smith   PetscFunctionBegin;
72a5ce6ee0Svictorle   if (lda>A->m) {
734ce68768SBarry Smith     nz = (PetscBLASInt)A->m;
74a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
75*f4df32b1SMatthew Knepley       BLASscal_(&nz,&oalpha,a->v+j*lda,&one);
76a5ce6ee0Svictorle     }
77a5ce6ee0Svictorle   } else {
784ce68768SBarry Smith     nz = (PetscBLASInt)A->m*A->n;
79*f4df32b1SMatthew Knepley     BLASscal_(&nz,&oalpha,a->v,&one);
80a5ce6ee0Svictorle   }
81efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
823a40ed3dSBarry Smith   PetscFunctionReturn(0);
8380cd9d93SLois Curfman McInnes }
8480cd9d93SLois Curfman McInnes 
85289bc588SBarry Smith /* ---------------------------------------------------------------*/
866831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
876831982aSBarry Smith    rather than put it in the Mat->row slot.*/
884a2ae208SSatish Balay #undef __FUNCT__
894a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
90dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
91289bc588SBarry Smith {
92a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
9381f479a6Svictorle   PetscFunctionBegin;
94a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
95a07ab388SBarry Smith #else
96c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
976849ba73SBarry Smith   PetscErrorCode ierr;
984ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,m = (PetscBLASInt)A->m,info;
9955659b69SBarry Smith 
1003a40ed3dSBarry Smith   PetscFunctionBegin;
101289bc588SBarry Smith   if (!mat->pivots) {
1024ce68768SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
10352e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,A->m*sizeof(PetscBLASInt));CHKERRQ(ierr);
104289bc588SBarry Smith   }
105f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
106273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
10771044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
10829bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10929bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
110efee365bSSatish Balay   ierr = PetscLogFlops((2*A->n*A->n*A->n)/3);CHKERRQ(ierr);
111a07ab388SBarry Smith #endif
1123a40ed3dSBarry Smith   PetscFunctionReturn(0);
113289bc588SBarry Smith }
1146ee01492SSatish Balay 
1154a2ae208SSatish Balay #undef __FUNCT__
1164a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
117dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11802cad45dSBarry Smith {
11902cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1206849ba73SBarry Smith   PetscErrorCode ierr;
12113f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
12202cad45dSBarry Smith   Mat            newi;
12302cad45dSBarry Smith 
1243a40ed3dSBarry Smith   PetscFunctionBegin;
125f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&newi);CHKERRQ(ierr);
126f69a0ea3SMatthew Knepley   ierr = MatSetSizes(newi,A->m,A->n,A->m,A->n);CHKERRQ(ierr);
1275c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1285c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1295609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
130a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
131a5ce6ee0Svictorle     if (lda>A->m) {
132a5ce6ee0Svictorle       m = A->m;
133a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
134a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
135a5ce6ee0Svictorle       }
136a5ce6ee0Svictorle     } else {
13787828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
13802cad45dSBarry Smith     }
139a5ce6ee0Svictorle   }
14002cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
14102cad45dSBarry Smith   *newmat = newi;
1423a40ed3dSBarry Smith   PetscFunctionReturn(0);
14302cad45dSBarry Smith }
14402cad45dSBarry Smith 
1454a2ae208SSatish Balay #undef __FUNCT__
1464a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
147dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
148289bc588SBarry Smith {
149dfbe8321SBarry Smith   PetscErrorCode ierr;
1503a40ed3dSBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
1525609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154289bc588SBarry Smith }
1556ee01492SSatish Balay 
1564a2ae208SSatish Balay #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
158af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
159289bc588SBarry Smith {
16002cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1616849ba73SBarry Smith   PetscErrorCode ierr;
16213f74950SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j;
1634482741eSBarry Smith   MatFactorInfo  info;
1643a40ed3dSBarry Smith 
1653a40ed3dSBarry Smith   PetscFunctionBegin;
16602cad45dSBarry Smith   /* copy the numerical values */
1671b807ce4Svictorle   if (lda1>m || lda2>m ) {
1681b807ce4Svictorle     for (j=0; j<n; j++) {
1691b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1701b807ce4Svictorle     }
1711b807ce4Svictorle   } else {
17287828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1731b807ce4Svictorle   }
17402cad45dSBarry Smith   (*fact)->factor = 0;
1754482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177289bc588SBarry Smith }
1786ee01492SSatish Balay 
1794a2ae208SSatish Balay #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
181dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
182289bc588SBarry Smith {
183dfbe8321SBarry Smith   PetscErrorCode ierr;
1843a40ed3dSBarry Smith 
1853a40ed3dSBarry Smith   PetscFunctionBegin;
186ceb03754SKris Buschelman   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr);
1873a40ed3dSBarry Smith   PetscFunctionReturn(0);
188289bc588SBarry Smith }
1896ee01492SSatish Balay 
1904a2ae208SSatish Balay #undef __FUNCT__
1914a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
192dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
193289bc588SBarry Smith {
194a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
195a07ab388SBarry Smith   PetscFunctionBegin;
196a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
197a07ab388SBarry Smith #else
198c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1996849ba73SBarry Smith   PetscErrorCode ierr;
2004ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,info;
20155659b69SBarry Smith 
2023a40ed3dSBarry Smith   PetscFunctionBegin;
203752f5567SLois Curfman McInnes   if (mat->pivots) {
204606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
20552e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-A->m*sizeof(PetscInt));CHKERRQ(ierr);
206752f5567SLois Curfman McInnes     mat->pivots = 0;
207752f5567SLois Curfman McInnes   }
208273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
20971044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
21077431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
211c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
212efee365bSSatish Balay   ierr = PetscLogFlops((A->n*A->n*A->n)/3);CHKERRQ(ierr);
213a07ab388SBarry Smith #endif
2143a40ed3dSBarry Smith   PetscFunctionReturn(0);
215289bc588SBarry Smith }
216289bc588SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2180b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
219af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2200b4b3355SBarry Smith {
221dfbe8321SBarry Smith   PetscErrorCode ierr;
22215e8a5b3SHong Zhang   MatFactorInfo  info;
2230b4b3355SBarry Smith 
2240b4b3355SBarry Smith   PetscFunctionBegin;
22515e8a5b3SHong Zhang   info.fill = 1.0;
22615e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2270b4b3355SBarry Smith   PetscFunctionReturn(0);
2280b4b3355SBarry Smith }
2290b4b3355SBarry Smith 
2300b4b3355SBarry Smith #undef __FUNCT__
2314a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
232dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
233289bc588SBarry Smith {
234c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2356849ba73SBarry Smith   PetscErrorCode ierr;
2364ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, one = 1,info;
23787828ca2SBarry Smith   PetscScalar    *x,*y;
23867e560aaSBarry Smith 
2393a40ed3dSBarry Smith   PetscFunctionBegin;
240273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2411ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2421ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
24387828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
244c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
245ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
24680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
247ae7cfcebSSatish Balay #else
24871044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2494ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
250ae7cfcebSSatish Balay #endif
2517a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
252ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
25380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
254ae7cfcebSSatish Balay #else
25571044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2564ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
257ae7cfcebSSatish Balay #endif
258289bc588SBarry Smith   }
25929bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2601ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2611ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
262efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2633a40ed3dSBarry Smith   PetscFunctionReturn(0);
264289bc588SBarry Smith }
2656ee01492SSatish Balay 
2664a2ae208SSatish Balay #undef __FUNCT__
2674a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
268dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
269da3a660dSBarry Smith {
270c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
271dfbe8321SBarry Smith   PetscErrorCode ierr;
2724ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->m,one = 1,info;
27387828ca2SBarry Smith   PetscScalar    *x,*y;
27467e560aaSBarry Smith 
2753a40ed3dSBarry Smith   PetscFunctionBegin;
276273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2771ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2781ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
27987828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
280752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
281da3a660dSBarry Smith   if (mat->pivots) {
282ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
28380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
284ae7cfcebSSatish Balay #else
28571044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2864ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
287ae7cfcebSSatish Balay #endif
2887a97a34bSBarry Smith   } else {
289ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
29080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
291ae7cfcebSSatish Balay #else
29271044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2934ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
294ae7cfcebSSatish Balay #endif
295da3a660dSBarry Smith   }
2961ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2971ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
298efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr);
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
300da3a660dSBarry Smith }
3016ee01492SSatish Balay 
3024a2ae208SSatish Balay #undef __FUNCT__
3034a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
304dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
305da3a660dSBarry Smith {
306c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
307dfbe8321SBarry Smith   PetscErrorCode ierr;
3084ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
30987828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
310da3a660dSBarry Smith   Vec            tmp = 0;
31167e560aaSBarry Smith 
3123a40ed3dSBarry Smith   PetscFunctionBegin;
3131ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3141ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
315273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
316da3a660dSBarry Smith   if (yy == zz) {
31778b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
31852e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
31978b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
320da3a660dSBarry Smith   }
32187828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
322752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
323da3a660dSBarry Smith   if (mat->pivots) {
324ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
32580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
326ae7cfcebSSatish Balay #else
32771044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3282ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
329ae7cfcebSSatish Balay #endif
330a8c6a408SBarry Smith   } else {
331ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
33280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
333ae7cfcebSSatish Balay #else
33471044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3352ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
336ae7cfcebSSatish Balay #endif
337da3a660dSBarry Smith   }
3382dcb1b2aSMatthew Knepley   if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
3392dcb1b2aSMatthew Knepley   else     {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);}
3401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3411ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
342efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
344da3a660dSBarry Smith }
34567e560aaSBarry Smith 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
348dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
349da3a660dSBarry Smith {
350c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3516849ba73SBarry Smith   PetscErrorCode ierr;
3524ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
35387828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
354da3a660dSBarry Smith   Vec            tmp;
35567e560aaSBarry Smith 
3563a40ed3dSBarry Smith   PetscFunctionBegin;
357273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3591ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
360da3a660dSBarry Smith   if (yy == zz) {
36178b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
36252e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
36378b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
364da3a660dSBarry Smith   }
36587828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
366752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
367da3a660dSBarry Smith   if (mat->pivots) {
368ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
370ae7cfcebSSatish Balay #else
37171044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3722ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
373ae7cfcebSSatish Balay #endif
3743a40ed3dSBarry Smith   } else {
375ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
37680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
377ae7cfcebSSatish Balay #else
37871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3792ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
380ae7cfcebSSatish Balay #endif
381da3a660dSBarry Smith   }
38290f02eecSBarry Smith   if (tmp) {
3832dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
38490f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3853a40ed3dSBarry Smith   } else {
3862dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
38790f02eecSBarry Smith   }
3881ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3891ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
390efee365bSSatish Balay   ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr);
3913a40ed3dSBarry Smith   PetscFunctionReturn(0);
392da3a660dSBarry Smith }
393289bc588SBarry Smith /* ------------------------------------------------------------------*/
3944a2ae208SSatish Balay #undef __FUNCT__
3954a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
39613f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
397289bc588SBarry Smith {
398c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
39987828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
400dfbe8321SBarry Smith   PetscErrorCode ierr;
40113f74950SBarry Smith   PetscInt       m = A->m,i;
402aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
4034ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
404bc1b551cSSatish Balay #endif
405289bc588SBarry Smith 
4063a40ed3dSBarry Smith   PetscFunctionBegin;
407289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
40871044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
4092dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
410289bc588SBarry Smith   }
4111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4121ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
413b965ef7fSBarry Smith   its  = its*lits;
41477431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
415289bc588SBarry Smith   while (its--) {
416fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
417289bc588SBarry Smith       for (i=0; i<m; i++) {
418aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
419f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
420f1747703SBarry Smith            not happy about returning a double complex */
42113f74950SBarry Smith         PetscInt         _i;
42287828ca2SBarry Smith         PetscScalar sum = b[i];
423f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4243f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
425f1747703SBarry Smith         }
426f1747703SBarry Smith         xt = sum;
427f1747703SBarry Smith #else
42871044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
429f1747703SBarry Smith #endif
43055a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
431289bc588SBarry Smith       }
432289bc588SBarry Smith     }
433fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
434289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
435aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
436f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
437f1747703SBarry Smith            not happy about returning a double complex */
43813f74950SBarry Smith         PetscInt         _i;
43987828ca2SBarry Smith         PetscScalar sum = b[i];
440f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4413f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
442f1747703SBarry Smith         }
443f1747703SBarry Smith         xt = sum;
444f1747703SBarry Smith #else
44571044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
446f1747703SBarry Smith #endif
44755a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
448289bc588SBarry Smith       }
449289bc588SBarry Smith     }
450289bc588SBarry Smith   }
4511ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4521ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4533a40ed3dSBarry Smith   PetscFunctionReturn(0);
454289bc588SBarry Smith }
455289bc588SBarry Smith 
456289bc588SBarry Smith /* -----------------------------------------------------------------*/
4574a2ae208SSatish Balay #undef __FUNCT__
4584a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
459dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
460289bc588SBarry Smith {
461c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
46287828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
463dfbe8321SBarry Smith   PetscErrorCode ierr;
4644ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n,_One=1;
465ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4663a40ed3dSBarry Smith 
4673a40ed3dSBarry Smith   PetscFunctionBegin;
468273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4691ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4701ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
47171044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4731ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
474efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->n);CHKERRQ(ierr);
4753a40ed3dSBarry Smith   PetscFunctionReturn(0);
476289bc588SBarry Smith }
4776ee01492SSatish Balay 
4784a2ae208SSatish Balay #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
480dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
481289bc588SBarry Smith {
482c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
48387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
484dfbe8321SBarry Smith   PetscErrorCode ierr;
4854ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
4863a40ed3dSBarry Smith 
4873a40ed3dSBarry Smith   PetscFunctionBegin;
488273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4901ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
49171044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4921ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4931ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
494efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n - A->m);CHKERRQ(ierr);
4953a40ed3dSBarry Smith   PetscFunctionReturn(0);
496289bc588SBarry Smith }
4976ee01492SSatish Balay 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
500dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
501289bc588SBarry Smith {
502c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
50387828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
504dfbe8321SBarry Smith   PetscErrorCode ierr;
5054ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
5063a40ed3dSBarry Smith 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
508273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
509600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5101ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5111ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
51271044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5131ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5141ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
515efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5163a40ed3dSBarry Smith   PetscFunctionReturn(0);
517289bc588SBarry Smith }
5186ee01492SSatish Balay 
5194a2ae208SSatish Balay #undef __FUNCT__
5204a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
521dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
522289bc588SBarry Smith {
523c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
52487828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
525dfbe8321SBarry Smith   PetscErrorCode ierr;
5264ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
52787828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5283a40ed3dSBarry Smith 
5293a40ed3dSBarry Smith   PetscFunctionBegin;
530273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
531600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5321ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5331ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
53471044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5361ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
537efee365bSSatish Balay   ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr);
5383a40ed3dSBarry Smith   PetscFunctionReturn(0);
539289bc588SBarry Smith }
540289bc588SBarry Smith 
541289bc588SBarry Smith /* -----------------------------------------------------------------*/
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
54413f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
545289bc588SBarry Smith {
546c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
54787828ca2SBarry Smith   PetscScalar    *v;
5486849ba73SBarry Smith   PetscErrorCode ierr;
54913f74950SBarry Smith   PetscInt       i;
55067e560aaSBarry Smith 
5513a40ed3dSBarry Smith   PetscFunctionBegin;
552273d9f13SBarry Smith   *ncols = A->n;
553289bc588SBarry Smith   if (cols) {
55413f74950SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
555273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
556289bc588SBarry Smith   }
557289bc588SBarry Smith   if (vals) {
55887828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
559289bc588SBarry Smith     v    = mat->v + row;
5601b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
561289bc588SBarry Smith   }
5623a40ed3dSBarry Smith   PetscFunctionReturn(0);
563289bc588SBarry Smith }
5646ee01492SSatish Balay 
5654a2ae208SSatish Balay #undef __FUNCT__
5664a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
56713f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
568289bc588SBarry Smith {
569dfbe8321SBarry Smith   PetscErrorCode ierr;
570606d414cSSatish Balay   PetscFunctionBegin;
571606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
572606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5733a40ed3dSBarry Smith   PetscFunctionReturn(0);
574289bc588SBarry Smith }
575289bc588SBarry Smith /* ----------------------------------------------------------------*/
5764a2ae208SSatish Balay #undef __FUNCT__
5774a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
57813f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
579289bc588SBarry Smith {
580c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
58113f74950SBarry Smith   PetscInt     i,j,idx=0;
582d6dfbf8fSBarry Smith 
5833a40ed3dSBarry Smith   PetscFunctionBegin;
584289bc588SBarry Smith   if (!mat->roworiented) {
585dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
586289bc588SBarry Smith       for (j=0; j<n; j++) {
587cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5882515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58977431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
59058804f6dSBarry Smith #endif
591289bc588SBarry Smith         for (i=0; i<m; i++) {
592cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5932515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59477431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
59558804f6dSBarry Smith #endif
596cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
597289bc588SBarry Smith         }
598289bc588SBarry Smith       }
5993a40ed3dSBarry Smith     } else {
600289bc588SBarry Smith       for (j=0; j<n; j++) {
601cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
6022515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60377431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
60458804f6dSBarry Smith #endif
605289bc588SBarry Smith         for (i=0; i<m; i++) {
606cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6072515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60877431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
60958804f6dSBarry Smith #endif
610cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
611289bc588SBarry Smith         }
612289bc588SBarry Smith       }
613289bc588SBarry Smith     }
6143a40ed3dSBarry Smith   } else {
615dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
616e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
617cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6182515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
61977431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
62058804f6dSBarry Smith #endif
621e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
622cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6232515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
62477431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
62558804f6dSBarry Smith #endif
626cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
627e8d4e0b9SBarry Smith         }
628e8d4e0b9SBarry Smith       }
6293a40ed3dSBarry Smith     } else {
630289bc588SBarry Smith       for (i=0; i<m; i++) {
631cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6322515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63377431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
63458804f6dSBarry Smith #endif
635289bc588SBarry Smith         for (j=0; j<n; j++) {
636cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6372515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63877431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
63958804f6dSBarry Smith #endif
640cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
641289bc588SBarry Smith         }
642289bc588SBarry Smith       }
643289bc588SBarry Smith     }
644e8d4e0b9SBarry Smith   }
6453a40ed3dSBarry Smith   PetscFunctionReturn(0);
646289bc588SBarry Smith }
647e8d4e0b9SBarry Smith 
6484a2ae208SSatish Balay #undef __FUNCT__
6494a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
65013f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
651ae80bb75SLois Curfman McInnes {
652ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
65313f74950SBarry Smith   PetscInt     i,j;
65487828ca2SBarry Smith   PetscScalar  *vpt = v;
655ae80bb75SLois Curfman McInnes 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
657ae80bb75SLois Curfman McInnes   /* row-oriented output */
658ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
659ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6601b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
661ae80bb75SLois Curfman McInnes     }
662ae80bb75SLois Curfman McInnes   }
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
664ae80bb75SLois Curfman McInnes }
665ae80bb75SLois Curfman McInnes 
666289bc588SBarry Smith /* -----------------------------------------------------------------*/
667289bc588SBarry Smith 
668e090d566SSatish Balay #include "petscsys.h"
66988e32edaSLois Curfman McInnes 
6704a2ae208SSatish Balay #undef __FUNCT__
6714a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
672f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A)
67388e32edaSLois Curfman McInnes {
67488e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
67588e32edaSLois Curfman McInnes   Mat            B;
6766849ba73SBarry Smith   PetscErrorCode ierr;
67713f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
67813f74950SBarry Smith   int            fd;
67913f74950SBarry Smith   PetscMPIInt    size;
68013f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
68187828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
68219bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
68388e32edaSLois Curfman McInnes 
6843a40ed3dSBarry Smith   PetscFunctionBegin;
685d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
68629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
687b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6880752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
689552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
69088e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
69188e32edaSLois Curfman McInnes 
69290ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
693f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
694f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
695be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6965c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
69790ace30eSBarry Smith     B    = *A;
69890ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6994a41a4d3SSatish Balay     v    = a->v;
7004a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
7014a41a4d3SSatish Balay        from row major to column major */
702b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
70390ace30eSBarry Smith     /* read in nonzero values */
7044a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
7054a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
7064a41a4d3SSatish Balay     for (j=0; j<N; j++) {
7074a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7084a41a4d3SSatish Balay         *v++ =w[i*N+j];
7094a41a4d3SSatish Balay       }
7104a41a4d3SSatish Balay     }
711606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7126d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7136d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71490ace30eSBarry Smith   } else {
71588e32edaSLois Curfman McInnes     /* read row lengths */
71613f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7170752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
71888e32edaSLois Curfman McInnes 
71988e32edaSLois Curfman McInnes     /* create our matrix */
720f69a0ea3SMatthew Knepley     ierr = MatCreate(comm,A);CHKERRQ(ierr);
721f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr);
722be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7235c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
72488e32edaSLois Curfman McInnes     B = *A;
72588e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
72688e32edaSLois Curfman McInnes     v = a->v;
72788e32edaSLois Curfman McInnes 
72888e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
72913f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
730b0a32e0cSBarry Smith     cols = scols;
7310752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
73287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
733b0a32e0cSBarry Smith     vals = svals;
7340752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
73588e32edaSLois Curfman McInnes 
73688e32edaSLois Curfman McInnes     /* insert into matrix */
73788e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
73888e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
73988e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
74088e32edaSLois Curfman McInnes     }
741606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
742606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
743606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
74488e32edaSLois Curfman McInnes 
7456d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7466d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74790ace30eSBarry Smith   }
7483a40ed3dSBarry Smith   PetscFunctionReturn(0);
74988e32edaSLois Curfman McInnes }
75088e32edaSLois Curfman McInnes 
751e090d566SSatish Balay #include "petscsys.h"
752932b0c3eSLois Curfman McInnes 
7534a2ae208SSatish Balay #undef __FUNCT__
7544a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7556849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
756289bc588SBarry Smith {
757932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
758dfbe8321SBarry Smith   PetscErrorCode    ierr;
75913f74950SBarry Smith   PetscInt          i,j;
7602dcb1b2aSMatthew Knepley   const char        *name;
76187828ca2SBarry Smith   PetscScalar       *v;
762f3ef73ceSBarry Smith   PetscViewerFormat format;
763932b0c3eSLois Curfman McInnes 
7643a40ed3dSBarry Smith   PetscFunctionBegin;
765435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
766b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
767456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7683a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
769fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
770b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
771273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
77244cd7ae7SLois Curfman McInnes       v = a->v + i;
77377431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
774273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
775aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
776329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
77777431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
778329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
77977431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7806831982aSBarry Smith         }
78180cd9d93SLois Curfman McInnes #else
7826831982aSBarry Smith         if (*v) {
78377431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,*v);CHKERRQ(ierr);
7846831982aSBarry Smith         }
78580cd9d93SLois Curfman McInnes #endif
7861b807ce4Svictorle         v += a->lda;
78780cd9d93SLois Curfman McInnes       }
788b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78980cd9d93SLois Curfman McInnes     }
790b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7913a40ed3dSBarry Smith   } else {
792b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
793aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
794ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
79547989497SBarry Smith     /* determine if matrix has all real values */
79647989497SBarry Smith     v = a->v;
797201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
798ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
79947989497SBarry Smith     }
80047989497SBarry Smith #endif
801fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
8023a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
80377431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->m,A->n);CHKERRQ(ierr);
80477431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->m,A->n);CHKERRQ(ierr);
805fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
806ffac6cdbSBarry Smith     }
807ffac6cdbSBarry Smith 
808273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
809932b0c3eSLois Curfman McInnes       v = a->v + i;
810273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
811aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
81247989497SBarry Smith         if (allreal) {
8131b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
81447989497SBarry Smith         } else {
8151b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
81647989497SBarry Smith         }
817289bc588SBarry Smith #else
8181b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
819289bc588SBarry Smith #endif
8201b807ce4Svictorle 	v += a->lda;
821289bc588SBarry Smith       }
822b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
823289bc588SBarry Smith     }
824fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
825b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
826ffac6cdbSBarry Smith     }
827b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
828da3a660dSBarry Smith   }
829b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8303a40ed3dSBarry Smith   PetscFunctionReturn(0);
831289bc588SBarry Smith }
832289bc588SBarry Smith 
8334a2ae208SSatish Balay #undef __FUNCT__
8344a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8356849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
836932b0c3eSLois Curfman McInnes {
837932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8386849ba73SBarry Smith   PetscErrorCode    ierr;
83913f74950SBarry Smith   int               fd;
84013f74950SBarry Smith   PetscInt          ict,j,n = A->n,m = A->m,i,*col_lens,nz = m*n;
84187828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
842f3ef73ceSBarry Smith   PetscViewerFormat format;
843932b0c3eSLois Curfman McInnes 
8443a40ed3dSBarry Smith   PetscFunctionBegin;
845b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
84690ace30eSBarry Smith 
847b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
848fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
84990ace30eSBarry Smith     /* store the matrix as a dense matrix */
85013f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
851552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
85290ace30eSBarry Smith     col_lens[1] = m;
85390ace30eSBarry Smith     col_lens[2] = n;
85490ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8556f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
856606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
85790ace30eSBarry Smith 
85890ace30eSBarry Smith     /* write out matrix, by rows */
85987828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
86090ace30eSBarry Smith     v    = a->v;
86190ace30eSBarry Smith     for (i=0; i<m; i++) {
86290ace30eSBarry Smith       for (j=0; j<n; j++) {
86390ace30eSBarry Smith         vals[i + j*m] = *v++;
86490ace30eSBarry Smith       }
86590ace30eSBarry Smith     }
8666f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
867606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
86890ace30eSBarry Smith   } else {
86913f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
870552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
871932b0c3eSLois Curfman McInnes     col_lens[1] = m;
872932b0c3eSLois Curfman McInnes     col_lens[2] = n;
873932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
874932b0c3eSLois Curfman McInnes 
875932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
876932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8776f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
878932b0c3eSLois Curfman McInnes 
879932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
880932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
881932b0c3eSLois Curfman McInnes     ict = 0;
882932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
883932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
884932b0c3eSLois Curfman McInnes     }
8856f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
886606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
887932b0c3eSLois Curfman McInnes 
888932b0c3eSLois Curfman McInnes     /* store nonzero values */
88987828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
890932b0c3eSLois Curfman McInnes     ict  = 0;
891932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
892932b0c3eSLois Curfman McInnes       v = a->v + i;
893932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8941b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
895932b0c3eSLois Curfman McInnes       }
896932b0c3eSLois Curfman McInnes     }
8976f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
898606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
89990ace30eSBarry Smith   }
9003a40ed3dSBarry Smith   PetscFunctionReturn(0);
901932b0c3eSLois Curfman McInnes }
902932b0c3eSLois Curfman McInnes 
9034a2ae208SSatish Balay #undef __FUNCT__
9044a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
905dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
906f1af5d2fSBarry Smith {
907f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
908f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9096849ba73SBarry Smith   PetscErrorCode    ierr;
91013f74950SBarry Smith   PetscInt          m = A->m,n = A->n,color,i,j;
91187828ca2SBarry Smith   PetscScalar       *v = a->v;
912b0a32e0cSBarry Smith   PetscViewer       viewer;
913b0a32e0cSBarry Smith   PetscDraw         popup;
914329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
915f3ef73ceSBarry Smith   PetscViewerFormat format;
916f1af5d2fSBarry Smith 
917f1af5d2fSBarry Smith   PetscFunctionBegin;
918f1af5d2fSBarry Smith 
919f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
920b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
921b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
922f1af5d2fSBarry Smith 
923f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
924fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
925f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
926b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
927f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
928f1af5d2fSBarry Smith       x_l = j;
929f1af5d2fSBarry Smith       x_r = x_l + 1.0;
930f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
931f1af5d2fSBarry Smith         y_l = m - i - 1.0;
932f1af5d2fSBarry Smith         y_r = y_l + 1.0;
933f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
934329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
935b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
936329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
937b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
938f1af5d2fSBarry Smith         } else {
939f1af5d2fSBarry Smith           continue;
940f1af5d2fSBarry Smith         }
941f1af5d2fSBarry Smith #else
942f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
943b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
944f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
945b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
946f1af5d2fSBarry Smith         } else {
947f1af5d2fSBarry Smith           continue;
948f1af5d2fSBarry Smith         }
949f1af5d2fSBarry Smith #endif
950b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
951f1af5d2fSBarry Smith       }
952f1af5d2fSBarry Smith     }
953f1af5d2fSBarry Smith   } else {
954f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
955f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
956f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
957f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
958f1af5d2fSBarry Smith     }
959b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
960b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
961b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
962f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
963f1af5d2fSBarry Smith       x_l = j;
964f1af5d2fSBarry Smith       x_r = x_l + 1.0;
965f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
966f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
967f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
968b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
969b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
970f1af5d2fSBarry Smith       }
971f1af5d2fSBarry Smith     }
972f1af5d2fSBarry Smith   }
973f1af5d2fSBarry Smith   PetscFunctionReturn(0);
974f1af5d2fSBarry Smith }
975f1af5d2fSBarry Smith 
9764a2ae208SSatish Balay #undef __FUNCT__
9774a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
978dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
979f1af5d2fSBarry Smith {
980b0a32e0cSBarry Smith   PetscDraw      draw;
981f1af5d2fSBarry Smith   PetscTruth     isnull;
982329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
983dfbe8321SBarry Smith   PetscErrorCode ierr;
984f1af5d2fSBarry Smith 
985f1af5d2fSBarry Smith   PetscFunctionBegin;
986b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
987b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
988abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
989f1af5d2fSBarry Smith 
990f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
991273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
992f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
993b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
994b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
995f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
996f1af5d2fSBarry Smith   PetscFunctionReturn(0);
997f1af5d2fSBarry Smith }
998f1af5d2fSBarry Smith 
9994a2ae208SSatish Balay #undef __FUNCT__
10004a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
1001dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1002932b0c3eSLois Curfman McInnes {
1003dfbe8321SBarry Smith   PetscErrorCode ierr;
100432077d6dSBarry Smith   PetscTruth     issocket,iascii,isbinary,isdraw;
1005932b0c3eSLois Curfman McInnes 
10063a40ed3dSBarry Smith   PetscFunctionBegin;
1007b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
100832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1009fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1010fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10110f5bd95cSBarry Smith 
1012c45a1595SBarry Smith   if (iascii) {
1013c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10144cf18111SSatish Balay #if defined(PETSC_USE_SOCKET_VIEWER)
1015c45a1595SBarry Smith   } else if (issocket) {
10164cf18111SSatish Balay     Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1017634064b4SBarry Smith     if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1018b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
1019c45a1595SBarry Smith #endif
10200f5bd95cSBarry Smith   } else if (isbinary) {
10213a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1022f1af5d2fSBarry Smith   } else if (isdraw) {
1023f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10245cd90555SBarry Smith   } else {
1025958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1026932b0c3eSLois Curfman McInnes   }
10273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1028932b0c3eSLois Curfman McInnes }
1029289bc588SBarry Smith 
10304a2ae208SSatish Balay #undef __FUNCT__
10314a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1032dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1033289bc588SBarry Smith {
1034ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1035dfbe8321SBarry Smith   PetscErrorCode ierr;
103690f02eecSBarry Smith 
10373a40ed3dSBarry Smith   PetscFunctionBegin;
1038aa482453SBarry Smith #if defined(PETSC_USE_LOG)
103977431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->m,mat->n);
1040a5a9c739SBarry Smith #endif
1041606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1042606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1043606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1044901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1046289bc588SBarry Smith }
1047289bc588SBarry Smith 
10484a2ae208SSatish Balay #undef __FUNCT__
10494a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1050dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1051289bc588SBarry Smith {
1052c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10536849ba73SBarry Smith   PetscErrorCode ierr;
105413f74950SBarry Smith   PetscInt       k,j,m,n,M;
105587828ca2SBarry Smith   PetscScalar    *v,tmp;
105648b35521SBarry Smith 
10573a40ed3dSBarry Smith   PetscFunctionBegin;
10581b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10597c922b88SBarry Smith   if (!matout) { /* in place transpose */
1060a5ce6ee0Svictorle     if (m != n) {
1061634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1062d64ed03dSBarry Smith     } else {
1063d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1064289bc588SBarry Smith         for (k=0; k<j; k++) {
10651b807ce4Svictorle           tmp = v[j + k*M];
10661b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10671b807ce4Svictorle           v[k + j*M] = tmp;
1068289bc588SBarry Smith         }
1069289bc588SBarry Smith       }
1070d64ed03dSBarry Smith     }
10713a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1072d3e5ee88SLois Curfman McInnes     Mat          tmat;
1073ec8511deSBarry Smith     Mat_SeqDense *tmatd;
107487828ca2SBarry Smith     PetscScalar  *v2;
1075ea709b57SSatish Balay 
1076f69a0ea3SMatthew Knepley     ierr  = MatCreate(A->comm,&tmat);CHKERRQ(ierr);
1077f69a0ea3SMatthew Knepley     ierr  = MatSetSizes(tmat,A->n,A->m,A->n,A->m);CHKERRQ(ierr);
10785c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10795c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1080ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10810de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1082d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10831b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1084d3e5ee88SLois Curfman McInnes     }
10856d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10866d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1087d3e5ee88SLois Curfman McInnes     *matout = tmat;
108848b35521SBarry Smith   }
10893a40ed3dSBarry Smith   PetscFunctionReturn(0);
1090289bc588SBarry Smith }
1091289bc588SBarry Smith 
10924a2ae208SSatish Balay #undef __FUNCT__
10934a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1094dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1095289bc588SBarry Smith {
1096c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1097c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
109813f74950SBarry Smith   PetscInt     i,j;
109987828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
11009ea5d5aeSSatish Balay 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
1102273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1103273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11041b807ce4Svictorle   for (i=0; i<A1->m; i++) {
11051b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
11061b807ce4Svictorle     for (j=0; j<A1->n; j++) {
11073a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
11081b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
11091b807ce4Svictorle     }
1110289bc588SBarry Smith   }
111177c4ece6SBarry Smith   *flg = PETSC_TRUE;
11123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1113289bc588SBarry Smith }
1114289bc588SBarry Smith 
11154a2ae208SSatish Balay #undef __FUNCT__
11164a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1117dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1118289bc588SBarry Smith {
1119c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1120dfbe8321SBarry Smith   PetscErrorCode ierr;
112113f74950SBarry Smith   PetscInt       i,n,len;
112287828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
112344cd7ae7SLois Curfman McInnes 
11243a40ed3dSBarry Smith   PetscFunctionBegin;
11252dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
11267a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11271ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1128273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1129273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
113044cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11311b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1132289bc588SBarry Smith   }
11331ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1135289bc588SBarry Smith }
1136289bc588SBarry Smith 
11374a2ae208SSatish Balay #undef __FUNCT__
11384a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1139dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1140289bc588SBarry Smith {
1141c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
114287828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1143dfbe8321SBarry Smith   PetscErrorCode ierr;
114413f74950SBarry Smith   PetscInt       i,j,m = A->m,n = A->n;
114555659b69SBarry Smith 
11463a40ed3dSBarry Smith   PetscFunctionBegin;
114728988994SBarry Smith   if (ll) {
11487a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11491ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1150273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1151da3a660dSBarry Smith     for (i=0; i<m; i++) {
1152da3a660dSBarry Smith       x = l[i];
1153da3a660dSBarry Smith       v = mat->v + i;
1154da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1155da3a660dSBarry Smith     }
11561ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1157efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1158da3a660dSBarry Smith   }
115928988994SBarry Smith   if (rr) {
11607a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11611ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1162273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1163da3a660dSBarry Smith     for (i=0; i<n; i++) {
1164da3a660dSBarry Smith       x = r[i];
1165da3a660dSBarry Smith       v = mat->v + i*m;
1166da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1167da3a660dSBarry Smith     }
11681ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1169efee365bSSatish Balay     ierr = PetscLogFlops(n*m);CHKERRQ(ierr);
1170da3a660dSBarry Smith   }
11713a40ed3dSBarry Smith   PetscFunctionReturn(0);
1172289bc588SBarry Smith }
1173289bc588SBarry Smith 
11744a2ae208SSatish Balay #undef __FUNCT__
11754a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1176dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1177289bc588SBarry Smith {
1178c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
117987828ca2SBarry Smith   PetscScalar  *v = mat->v;
1180329f5518SBarry Smith   PetscReal    sum = 0.0;
118113f74950SBarry Smith   PetscInt     lda=mat->lda,m=A->m,i,j;
1182efee365bSSatish Balay   PetscErrorCode ierr;
118355659b69SBarry Smith 
11843a40ed3dSBarry Smith   PetscFunctionBegin;
1185289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1186a5ce6ee0Svictorle     if (lda>m) {
1187a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1188a5ce6ee0Svictorle 	v = mat->v+j*lda;
1189a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1190a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1191a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1192a5ce6ee0Svictorle #else
1193a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1194a5ce6ee0Svictorle #endif
1195a5ce6ee0Svictorle 	}
1196a5ce6ee0Svictorle       }
1197a5ce6ee0Svictorle     } else {
1198273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1199aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1200329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1201289bc588SBarry Smith #else
1202289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1203289bc588SBarry Smith #endif
1204289bc588SBarry Smith       }
1205a5ce6ee0Svictorle     }
1206064f8208SBarry Smith     *nrm = sqrt(sum);
1207efee365bSSatish Balay     ierr = PetscLogFlops(2*A->n*A->m);CHKERRQ(ierr);
12083a40ed3dSBarry Smith   } else if (type == NORM_1) {
1209064f8208SBarry Smith     *nrm = 0.0;
1210273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
12111b807ce4Svictorle       v = mat->v + j*mat->lda;
1212289bc588SBarry Smith       sum = 0.0;
1213273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
121433a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1215289bc588SBarry Smith       }
1216064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1217289bc588SBarry Smith     }
1218efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12193a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1220064f8208SBarry Smith     *nrm = 0.0;
1221273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1222289bc588SBarry Smith       v = mat->v + j;
1223289bc588SBarry Smith       sum = 0.0;
1224273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
12251b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1226289bc588SBarry Smith       }
1227064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1228289bc588SBarry Smith     }
1229efee365bSSatish Balay     ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr);
12303a40ed3dSBarry Smith   } else {
123129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1232289bc588SBarry Smith   }
12333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1234289bc588SBarry Smith }
1235289bc588SBarry Smith 
12364a2ae208SSatish Balay #undef __FUNCT__
12374a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1238dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1239289bc588SBarry Smith {
1240c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
124163ba0a88SBarry Smith   PetscErrorCode ierr;
124267e560aaSBarry Smith 
12433a40ed3dSBarry Smith   PetscFunctionBegin;
1244b5a2b587SKris Buschelman   switch (op) {
1245b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1246b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1247b5a2b587SKris Buschelman     break;
1248b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1249b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1250b5a2b587SKris Buschelman     break;
1251b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1252b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1253b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1254b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1255b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1256b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1257b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1258b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1259b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1260b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1261b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
126263ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_SeqDense:Option ignored\n"));CHKERRQ(ierr);
1263b5a2b587SKris Buschelman     break;
126477e54ba9SKris Buschelman   case MAT_SYMMETRIC:
126577e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12669a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12679a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12689a4540c5SBarry Smith   case MAT_HERMITIAN:
12699a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12709a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12719a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
127277e54ba9SKris Buschelman     break;
1273b5a2b587SKris Buschelman   default:
127429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12753a40ed3dSBarry Smith   }
12763a40ed3dSBarry Smith   PetscFunctionReturn(0);
1277289bc588SBarry Smith }
1278289bc588SBarry Smith 
12794a2ae208SSatish Balay #undef __FUNCT__
12804a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1281dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12826f0a148fSBarry Smith {
1283ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12846849ba73SBarry Smith   PetscErrorCode ierr;
128513f74950SBarry Smith   PetscInt       lda=l->lda,m=A->m,j;
12863a40ed3dSBarry Smith 
12873a40ed3dSBarry Smith   PetscFunctionBegin;
1288a5ce6ee0Svictorle   if (lda>m) {
1289a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1290a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1291a5ce6ee0Svictorle     }
1292a5ce6ee0Svictorle   } else {
129387828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1294a5ce6ee0Svictorle   }
12953a40ed3dSBarry Smith   PetscFunctionReturn(0);
12966f0a148fSBarry Smith }
12976f0a148fSBarry Smith 
12984a2ae208SSatish Balay #undef __FUNCT__
12994a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1300*f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
13016f0a148fSBarry Smith {
1302ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
1303*f4df32b1SMatthew Knepley   PetscInt       n = A->n,i,j;
130487828ca2SBarry Smith   PetscScalar    *slot;
130555659b69SBarry Smith 
13063a40ed3dSBarry Smith   PetscFunctionBegin;
13076f0a148fSBarry Smith   for (i=0; i<N; i++) {
13086f0a148fSBarry Smith     slot = l->v + rows[i];
13096f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
13106f0a148fSBarry Smith   }
1311*f4df32b1SMatthew Knepley   if (diag != 0.0) {
13126f0a148fSBarry Smith     for (i=0; i<N; i++) {
13136f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1314*f4df32b1SMatthew Knepley       *slot = diag;
13156f0a148fSBarry Smith     }
13166f0a148fSBarry Smith   }
13173a40ed3dSBarry Smith   PetscFunctionReturn(0);
13186f0a148fSBarry Smith }
1319557bce09SLois Curfman McInnes 
13204a2ae208SSatish Balay #undef __FUNCT__
13214a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1322dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
132364e87e97SBarry Smith {
1324c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13253a40ed3dSBarry Smith 
13263a40ed3dSBarry Smith   PetscFunctionBegin;
132764e87e97SBarry Smith   *array = mat->v;
13283a40ed3dSBarry Smith   PetscFunctionReturn(0);
132964e87e97SBarry Smith }
13300754003eSLois Curfman McInnes 
13314a2ae208SSatish Balay #undef __FUNCT__
13324a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1333dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1334ff14e315SSatish Balay {
13353a40ed3dSBarry Smith   PetscFunctionBegin;
133609b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13373a40ed3dSBarry Smith   PetscFunctionReturn(0);
1338ff14e315SSatish Balay }
13390754003eSLois Curfman McInnes 
13404a2ae208SSatish Balay #undef __FUNCT__
13414a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
134213f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13430754003eSLois Curfman McInnes {
1344c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13456849ba73SBarry Smith   PetscErrorCode ierr;
134613f74950SBarry Smith   PetscInt       i,j,m = A->m,*irow,*icol,nrows,ncols;
134787828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13480754003eSLois Curfman McInnes   Mat            newmat;
13490754003eSLois Curfman McInnes 
13503a40ed3dSBarry Smith   PetscFunctionBegin;
135178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
135278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1353e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1354e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13550754003eSLois Curfman McInnes 
1356182d2002SSatish Balay   /* Check submatrixcall */
1357182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
135813f74950SBarry Smith     PetscInt n_cols,n_rows;
1359182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
136029bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1361182d2002SSatish Balay     newmat = *B;
1362182d2002SSatish Balay   } else {
13630754003eSLois Curfman McInnes     /* Create and fill new matrix */
1364f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr);
1365f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
13665c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13675c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1368182d2002SSatish Balay   }
1369182d2002SSatish Balay 
1370182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1371182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1372182d2002SSatish Balay 
1373182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1374182d2002SSatish Balay     av = v + m*icol[i];
1375182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1376182d2002SSatish Balay       *bv++ = av[irow[j]];
13770754003eSLois Curfman McInnes     }
13780754003eSLois Curfman McInnes   }
1379182d2002SSatish Balay 
1380182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13816d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13826d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13830754003eSLois Curfman McInnes 
13840754003eSLois Curfman McInnes   /* Free work space */
138578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
138678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1387182d2002SSatish Balay   *B = newmat;
13883a40ed3dSBarry Smith   PetscFunctionReturn(0);
13890754003eSLois Curfman McInnes }
13900754003eSLois Curfman McInnes 
13914a2ae208SSatish Balay #undef __FUNCT__
13924a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
139313f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1394905e6a2fSBarry Smith {
13956849ba73SBarry Smith   PetscErrorCode ierr;
139613f74950SBarry Smith   PetscInt       i;
1397905e6a2fSBarry Smith 
13983a40ed3dSBarry Smith   PetscFunctionBegin;
1399905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1400b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1401905e6a2fSBarry Smith   }
1402905e6a2fSBarry Smith 
1403905e6a2fSBarry Smith   for (i=0; i<n; i++) {
14046a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1405905e6a2fSBarry Smith   }
14063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1407905e6a2fSBarry Smith }
1408905e6a2fSBarry Smith 
14094a2ae208SSatish Balay #undef __FUNCT__
14104a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1411dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14124b0e389bSBarry Smith {
14134b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14146849ba73SBarry Smith   PetscErrorCode ierr;
141513f74950SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j;
14163a40ed3dSBarry Smith 
14173a40ed3dSBarry Smith   PetscFunctionBegin;
141833f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
141933f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1420cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14213a40ed3dSBarry Smith     PetscFunctionReturn(0);
14223a40ed3dSBarry Smith   }
14230dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1424a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14250dbb7854Svictorle     for (j=0; j<n; j++) {
1426a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1427a5ce6ee0Svictorle     }
1428a5ce6ee0Svictorle   } else {
142987828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1430a5ce6ee0Svictorle   }
1431273d9f13SBarry Smith   PetscFunctionReturn(0);
1432273d9f13SBarry Smith }
1433273d9f13SBarry Smith 
14344a2ae208SSatish Balay #undef __FUNCT__
14354a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1436dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1437273d9f13SBarry Smith {
1438dfbe8321SBarry Smith   PetscErrorCode ierr;
1439273d9f13SBarry Smith 
1440273d9f13SBarry Smith   PetscFunctionBegin;
1441273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14423a40ed3dSBarry Smith   PetscFunctionReturn(0);
14434b0e389bSBarry Smith }
14444b0e389bSBarry Smith 
1445289bc588SBarry Smith /* -------------------------------------------------------------------*/
1446a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1447905e6a2fSBarry Smith        MatGetRow_SeqDense,
1448905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1449905e6a2fSBarry Smith        MatMult_SeqDense,
145097304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14517c922b88SBarry Smith        MatMultTranspose_SeqDense,
14527c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1453905e6a2fSBarry Smith        MatSolve_SeqDense,
1454905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14557c922b88SBarry Smith        MatSolveTranspose_SeqDense,
145697304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1457905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1458905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1459ec8511deSBarry Smith        MatRelax_SeqDense,
1460ec8511deSBarry Smith        MatTranspose_SeqDense,
146197304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1462905e6a2fSBarry Smith        MatEqual_SeqDense,
1463905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1464905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1465905e6a2fSBarry Smith        MatNorm_SeqDense,
146697304618SKris Buschelman /*20*/ 0,
1467905e6a2fSBarry Smith        0,
1468905e6a2fSBarry Smith        0,
1469905e6a2fSBarry Smith        MatSetOption_SeqDense,
1470905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
147197304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1472905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1473905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1474905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1475905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
147697304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1477273d9f13SBarry Smith        0,
1478905e6a2fSBarry Smith        0,
1479905e6a2fSBarry Smith        MatGetArray_SeqDense,
1480905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
148197304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1482a5ae1ecdSBarry Smith        0,
1483a5ae1ecdSBarry Smith        0,
1484a5ae1ecdSBarry Smith        0,
1485a5ae1ecdSBarry Smith        0,
148697304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1487a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1488a5ae1ecdSBarry Smith        0,
14894b0e389bSBarry Smith        MatGetValues_SeqDense,
1490a5ae1ecdSBarry Smith        MatCopy_SeqDense,
149197304618SKris Buschelman /*45*/ 0,
1492a5ae1ecdSBarry Smith        MatScale_SeqDense,
1493a5ae1ecdSBarry Smith        0,
1494a5ae1ecdSBarry Smith        0,
1495a5ae1ecdSBarry Smith        0,
1496521d7252SBarry Smith /*50*/ 0,
1497a5ae1ecdSBarry Smith        0,
1498a5ae1ecdSBarry Smith        0,
1499a5ae1ecdSBarry Smith        0,
1500a5ae1ecdSBarry Smith        0,
150197304618SKris Buschelman /*55*/ 0,
1502a5ae1ecdSBarry Smith        0,
1503a5ae1ecdSBarry Smith        0,
1504a5ae1ecdSBarry Smith        0,
1505a5ae1ecdSBarry Smith        0,
150697304618SKris Buschelman /*60*/ 0,
1507e03a110bSBarry Smith        MatDestroy_SeqDense,
1508e03a110bSBarry Smith        MatView_SeqDense,
150997304618SKris Buschelman        MatGetPetscMaps_Petsc,
151097304618SKris Buschelman        0,
151197304618SKris Buschelman /*65*/ 0,
151297304618SKris Buschelman        0,
151397304618SKris Buschelman        0,
151497304618SKris Buschelman        0,
151597304618SKris Buschelman        0,
151697304618SKris Buschelman /*70*/ 0,
151797304618SKris Buschelman        0,
151897304618SKris Buschelman        0,
151997304618SKris Buschelman        0,
152097304618SKris Buschelman        0,
152197304618SKris Buschelman /*75*/ 0,
152297304618SKris Buschelman        0,
152397304618SKris Buschelman        0,
152497304618SKris Buschelman        0,
152597304618SKris Buschelman        0,
152697304618SKris Buschelman /*80*/ 0,
152797304618SKris Buschelman        0,
152897304618SKris Buschelman        0,
152997304618SKris Buschelman        0,
1530865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1531865e5f61SKris Buschelman        0,
1532865e5f61SKris Buschelman        0,
1533865e5f61SKris Buschelman        0,
1534865e5f61SKris Buschelman        0,
1535865e5f61SKris Buschelman        0,
1536865e5f61SKris Buschelman /*90*/ 0,
1537865e5f61SKris Buschelman        0,
1538865e5f61SKris Buschelman        0,
1539865e5f61SKris Buschelman        0,
1540865e5f61SKris Buschelman        0,
1541865e5f61SKris Buschelman /*95*/ 0,
1542865e5f61SKris Buschelman        0,
1543865e5f61SKris Buschelman        0,
1544865e5f61SKris Buschelman        0};
154590ace30eSBarry Smith 
15464a2ae208SSatish Balay #undef __FUNCT__
15474a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
15484b828684SBarry Smith /*@C
1549fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1550d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1551d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1552289bc588SBarry Smith 
1553db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1554db81eaa0SLois Curfman McInnes 
155520563c6bSBarry Smith    Input Parameters:
1556db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15570c775827SLois Curfman McInnes .  m - number of rows
155818f449edSLois Curfman McInnes .  n - number of columns
1559db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1560dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
156120563c6bSBarry Smith 
156220563c6bSBarry Smith    Output Parameter:
156344cd7ae7SLois Curfman McInnes .  A - the matrix
156420563c6bSBarry Smith 
1565b259b22eSLois Curfman McInnes    Notes:
156618f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
156718f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1568b4fd4287SBarry Smith    set data=PETSC_NULL.
156918f449edSLois Curfman McInnes 
1570027ccd11SLois Curfman McInnes    Level: intermediate
1571027ccd11SLois Curfman McInnes 
1572dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1573d65003e9SLois Curfman McInnes 
1574db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
157520563c6bSBarry Smith @*/
1576be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1577289bc588SBarry Smith {
1578dfbe8321SBarry Smith   PetscErrorCode ierr;
15793b2fbd54SBarry Smith 
15803a40ed3dSBarry Smith   PetscFunctionBegin;
1581f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1582f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1583273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1584273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1585273d9f13SBarry Smith   PetscFunctionReturn(0);
1586273d9f13SBarry Smith }
1587273d9f13SBarry Smith 
15884a2ae208SSatish Balay #undef __FUNCT__
15894a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1590273d9f13SBarry Smith /*@C
1591273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1592273d9f13SBarry Smith 
1593273d9f13SBarry Smith    Collective on MPI_Comm
1594273d9f13SBarry Smith 
1595273d9f13SBarry Smith    Input Parameters:
1596273d9f13SBarry Smith +  A - the matrix
1597273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1598273d9f13SBarry Smith 
1599273d9f13SBarry Smith    Notes:
1600273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1601273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1602273d9f13SBarry Smith    set data=PETSC_NULL.
1603273d9f13SBarry Smith 
1604273d9f13SBarry Smith    Level: intermediate
1605273d9f13SBarry Smith 
1606273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1607273d9f13SBarry Smith 
1608273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1609273d9f13SBarry Smith @*/
1610be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1611273d9f13SBarry Smith {
1612dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1613a23d5eceSKris Buschelman 
1614a23d5eceSKris Buschelman   PetscFunctionBegin;
1615a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1616a23d5eceSKris Buschelman   if (f) {
1617a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1618a23d5eceSKris Buschelman   }
1619a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1620a23d5eceSKris Buschelman }
1621a23d5eceSKris Buschelman 
1622a23d5eceSKris Buschelman EXTERN_C_BEGIN
1623a23d5eceSKris Buschelman #undef __FUNCT__
1624a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1625be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1626a23d5eceSKris Buschelman {
1627273d9f13SBarry Smith   Mat_SeqDense   *b;
1628dfbe8321SBarry Smith   PetscErrorCode ierr;
1629273d9f13SBarry Smith 
1630273d9f13SBarry Smith   PetscFunctionBegin;
1631273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1632273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1633273d9f13SBarry Smith   if (!data) {
163487828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
163587828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1636273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
163752e6d16bSBarry Smith     ierr = PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));CHKERRQ(ierr);
1638273d9f13SBarry Smith   } else { /* user-allocated storage */
1639273d9f13SBarry Smith     b->v          = data;
1640273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1641273d9f13SBarry Smith   }
1642273d9f13SBarry Smith   PetscFunctionReturn(0);
1643273d9f13SBarry Smith }
1644a23d5eceSKris Buschelman EXTERN_C_END
1645273d9f13SBarry Smith 
16461b807ce4Svictorle #undef __FUNCT__
16471b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
16481b807ce4Svictorle /*@C
16491b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
16501b807ce4Svictorle 
16511b807ce4Svictorle   Input parameter:
16521b807ce4Svictorle + A - the matrix
16531b807ce4Svictorle - lda - the leading dimension
16541b807ce4Svictorle 
16551b807ce4Svictorle   Notes:
16561b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16571b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16581b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16591b807ce4Svictorle 
16601b807ce4Svictorle   Level: intermediate
16611b807ce4Svictorle 
16621b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16631b807ce4Svictorle 
16641b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16651b807ce4Svictorle @*/
1666be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda)
16671b807ce4Svictorle {
16681b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16691b807ce4Svictorle   PetscFunctionBegin;
167077431f27SBarry Smith   if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->m);
16711b807ce4Svictorle   b->lda = lda;
16721b807ce4Svictorle   PetscFunctionReturn(0);
16731b807ce4Svictorle }
16741b807ce4Svictorle 
16750bad9183SKris Buschelman /*MC
1676fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16770bad9183SKris Buschelman 
16780bad9183SKris Buschelman    Options Database Keys:
16790bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16800bad9183SKris Buschelman 
16810bad9183SKris Buschelman   Level: beginner
16820bad9183SKris Buschelman 
16830bad9183SKris Buschelman .seealso: MatCreateSeqDense
16840bad9183SKris Buschelman M*/
16850bad9183SKris Buschelman 
1686273d9f13SBarry Smith EXTERN_C_BEGIN
16874a2ae208SSatish Balay #undef __FUNCT__
16884a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1689be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B)
1690273d9f13SBarry Smith {
1691273d9f13SBarry Smith   Mat_SeqDense   *b;
1692dfbe8321SBarry Smith   PetscErrorCode ierr;
16937c334f02SBarry Smith   PetscMPIInt    size;
1694273d9f13SBarry Smith 
1695273d9f13SBarry Smith   PetscFunctionBegin;
1696273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
169729bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
169855659b69SBarry Smith 
1699273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1700273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1701273d9f13SBarry Smith 
1702b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1703549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
170444cd7ae7SLois Curfman McInnes   B->factor       = 0;
170590f02eecSBarry Smith   B->mapping      = 0;
170652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
170744cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
170818f449edSLois Curfman McInnes 
17098a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
17108a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1711a5ae1ecdSBarry Smith 
171244cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1713273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1714273d9f13SBarry Smith   b->v            = 0;
17151b807ce4Svictorle   b->lda          = B->m;
17164e220ebcSLois Curfman McInnes 
1717a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1718a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1719a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
17203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1721289bc588SBarry Smith }
1722273d9f13SBarry Smith EXTERN_C_END
1723