xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 52e6d16ba989af9362d0fcebb12e73dd7c9666ed)
167e560aaSBarry Smith /*
267e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
367e560aaSBarry Smith */
4289bc588SBarry Smith 
570f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
6b4c862e0SSatish Balay #include "petscblaslapack.h"
7289bc588SBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense"
10dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str)
111987afe7SBarry Smith {
121987afe7SBarry Smith   Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
1313f74950SBarry Smith   PetscInt     j;
144ce68768SBarry Smith   PetscBLASInt N = (PetscBLASInt)X->m*X->n,m=(PetscBLASInt)X->m,ldax = x->lda,lday=y->lda,one = 1;
153a40ed3dSBarry Smith 
163a40ed3dSBarry Smith   PetscFunctionBegin;
17a5ce6ee0Svictorle   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
18a5ce6ee0Svictorle   if (ldax>m || lday>m) {
19a5ce6ee0Svictorle     for (j=0; j<X->n; j++) {
2071044d3cSBarry Smith       BLASaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one);
21a5ce6ee0Svictorle     }
22a5ce6ee0Svictorle   } else {
2371044d3cSBarry Smith     BLASaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one);
24a5ce6ee0Svictorle   }
25b0a32e0cSBarry Smith   PetscLogFlops(2*N-1);
263a40ed3dSBarry Smith   PetscFunctionReturn(0);
271987afe7SBarry Smith }
281987afe7SBarry Smith 
294a2ae208SSatish Balay #undef __FUNCT__
304a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
31dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
32289bc588SBarry Smith {
334e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
3413f74950SBarry Smith   PetscInt     i,N = A->m*A->n,count = 0;
3587828ca2SBarry Smith   PetscScalar  *v = a->v;
363a40ed3dSBarry Smith 
373a40ed3dSBarry Smith   PetscFunctionBegin;
38289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
394e220ebcSLois Curfman McInnes 
40273d9f13SBarry Smith   info->rows_global       = (double)A->m;
41273d9f13SBarry Smith   info->columns_global    = (double)A->n;
42273d9f13SBarry Smith   info->rows_local        = (double)A->m;
43273d9f13SBarry Smith   info->columns_local     = (double)A->n;
444e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
454e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
464e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
474e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
484e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
494e220ebcSLois Curfman McInnes   info->mallocs           = 0;
504e220ebcSLois Curfman McInnes   info->memory            = A->mem;
514e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
524e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
534e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
544e220ebcSLois Curfman McInnes 
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56289bc588SBarry Smith }
57289bc588SBarry Smith 
584a2ae208SSatish Balay #undef __FUNCT__
594a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
60dfbe8321SBarry Smith PetscErrorCode MatScale_SeqDense(const PetscScalar *alpha,Mat A)
6180cd9d93SLois Curfman McInnes {
62273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
634ce68768SBarry Smith   PetscBLASInt one = 1,lda = a->lda,j,nz;
6480cd9d93SLois Curfman McInnes 
653a40ed3dSBarry Smith   PetscFunctionBegin;
66a5ce6ee0Svictorle   if (lda>A->m) {
674ce68768SBarry Smith     nz = (PetscBLASInt)A->m;
68a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
6971044d3cSBarry Smith       BLASscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one);
70a5ce6ee0Svictorle     }
71a5ce6ee0Svictorle   } else {
724ce68768SBarry Smith     nz = (PetscBLASInt)A->m*A->n;
7371044d3cSBarry Smith     BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one);
74a5ce6ee0Svictorle   }
75b0a32e0cSBarry Smith   PetscLogFlops(nz);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
7780cd9d93SLois Curfman McInnes }
7880cd9d93SLois Curfman McInnes 
79289bc588SBarry Smith /* ---------------------------------------------------------------*/
806831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
816831982aSBarry Smith    rather than put it in the Mat->row slot.*/
824a2ae208SSatish Balay #undef __FUNCT__
834a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
84dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
85289bc588SBarry Smith {
86a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
8781f479a6Svictorle   PetscFunctionBegin;
88a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
89a07ab388SBarry Smith #else
90c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
916849ba73SBarry Smith   PetscErrorCode ierr;
924ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,m = (PetscBLASInt)A->m,info;
9355659b69SBarry Smith 
943a40ed3dSBarry Smith   PetscFunctionBegin;
95289bc588SBarry Smith   if (!mat->pivots) {
964ce68768SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr);
97*52e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,A->m*sizeof(PetscBLASInt));CHKERRQ(ierr);
98289bc588SBarry Smith   }
99f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
100273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
10171044d3cSBarry Smith   LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info);
10229bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10329bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
104b0a32e0cSBarry Smith   PetscLogFlops((2*A->n*A->n*A->n)/3);
105a07ab388SBarry Smith #endif
1063a40ed3dSBarry Smith   PetscFunctionReturn(0);
107289bc588SBarry Smith }
1086ee01492SSatish Balay 
1094a2ae208SSatish Balay #undef __FUNCT__
1104a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
111dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11202cad45dSBarry Smith {
11302cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
1146849ba73SBarry Smith   PetscErrorCode ierr;
11513f74950SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
11602cad45dSBarry Smith   Mat            newi;
11702cad45dSBarry Smith 
1183a40ed3dSBarry Smith   PetscFunctionBegin;
1195c5985e7SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&newi);CHKERRQ(ierr);
1205c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1215c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1225609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
123a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
124a5ce6ee0Svictorle     if (lda>A->m) {
125a5ce6ee0Svictorle       m = A->m;
126a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
127a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
128a5ce6ee0Svictorle       }
129a5ce6ee0Svictorle     } else {
13087828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
13102cad45dSBarry Smith     }
132a5ce6ee0Svictorle   }
13302cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13402cad45dSBarry Smith   *newmat = newi;
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
13602cad45dSBarry Smith }
13702cad45dSBarry Smith 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
140dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
141289bc588SBarry Smith {
142dfbe8321SBarry Smith   PetscErrorCode ierr;
1433a40ed3dSBarry Smith 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
1455609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
147289bc588SBarry Smith }
1486ee01492SSatish Balay 
1494a2ae208SSatish Balay #undef __FUNCT__
1504a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
151af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
152289bc588SBarry Smith {
15302cad45dSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1546849ba73SBarry Smith   PetscErrorCode ierr;
15513f74950SBarry Smith   PetscInt       lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j;
1564482741eSBarry Smith   MatFactorInfo  info;
1573a40ed3dSBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
15902cad45dSBarry Smith   /* copy the numerical values */
1601b807ce4Svictorle   if (lda1>m || lda2>m ) {
1611b807ce4Svictorle     for (j=0; j<n; j++) {
1621b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1631b807ce4Svictorle     }
1641b807ce4Svictorle   } else {
16587828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1661b807ce4Svictorle   }
16702cad45dSBarry Smith   (*fact)->factor = 0;
1684482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1693a40ed3dSBarry Smith   PetscFunctionReturn(0);
170289bc588SBarry Smith }
1716ee01492SSatish Balay 
1724a2ae208SSatish Balay #undef __FUNCT__
1734a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
174dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
175289bc588SBarry Smith {
176dfbe8321SBarry Smith   PetscErrorCode ierr;
1773a40ed3dSBarry Smith 
1783a40ed3dSBarry Smith   PetscFunctionBegin;
1793a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
181289bc588SBarry Smith }
1826ee01492SSatish Balay 
1834a2ae208SSatish Balay #undef __FUNCT__
1844a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
185dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
186289bc588SBarry Smith {
187a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
188a07ab388SBarry Smith   PetscFunctionBegin;
189a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
190a07ab388SBarry Smith #else
191c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1926849ba73SBarry Smith   PetscErrorCode ierr;
1934ce68768SBarry Smith   PetscBLASInt   n = (PetscBLASInt)A->n,info;
19455659b69SBarry Smith 
1953a40ed3dSBarry Smith   PetscFunctionBegin;
196752f5567SLois Curfman McInnes   if (mat->pivots) {
197606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
198*52e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-A->m*sizeof(PetscInt));CHKERRQ(ierr);
199752f5567SLois Curfman McInnes     mat->pivots = 0;
200752f5567SLois Curfman McInnes   }
201273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
20271044d3cSBarry Smith   LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info);
20377431f27SBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
204c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
205b0a32e0cSBarry Smith   PetscLogFlops((A->n*A->n*A->n)/3);
206a07ab388SBarry Smith #endif
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
208289bc588SBarry Smith }
209289bc588SBarry Smith 
2104a2ae208SSatish Balay #undef __FUNCT__
2110b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
212af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact)
2130b4b3355SBarry Smith {
214dfbe8321SBarry Smith   PetscErrorCode ierr;
21515e8a5b3SHong Zhang   MatFactorInfo  info;
2160b4b3355SBarry Smith 
2170b4b3355SBarry Smith   PetscFunctionBegin;
21815e8a5b3SHong Zhang   info.fill = 1.0;
21915e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2200b4b3355SBarry Smith   PetscFunctionReturn(0);
2210b4b3355SBarry Smith }
2220b4b3355SBarry Smith 
2230b4b3355SBarry Smith #undef __FUNCT__
2244a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
225dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
226289bc588SBarry Smith {
227c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
2286849ba73SBarry Smith   PetscErrorCode ierr;
2294ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, one = 1,info;
23087828ca2SBarry Smith   PetscScalar    *x,*y;
23167e560aaSBarry Smith 
2323a40ed3dSBarry Smith   PetscFunctionBegin;
233273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2341ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2351ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
23687828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
237c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
238ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
23980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
240ae7cfcebSSatish Balay #else
24171044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2424ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve");
243ae7cfcebSSatish Balay #endif
2447a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
245ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
24680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
247ae7cfcebSSatish Balay #else
24871044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2494ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve");
250ae7cfcebSSatish Balay #endif
251289bc588SBarry Smith   }
25229bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2531ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2541ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
255b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
257289bc588SBarry Smith }
2586ee01492SSatish Balay 
2594a2ae208SSatish Balay #undef __FUNCT__
2604a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
261dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
262da3a660dSBarry Smith {
263c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
264dfbe8321SBarry Smith   PetscErrorCode ierr;
2654ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt) A->m,one = 1,info;
26687828ca2SBarry Smith   PetscScalar    *x,*y;
26767e560aaSBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionBegin;
269273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2701ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2711ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
27287828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
273752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
274da3a660dSBarry Smith   if (mat->pivots) {
275ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
27680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
277ae7cfcebSSatish Balay #else
27871044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
2794ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
280ae7cfcebSSatish Balay #endif
2817a97a34bSBarry Smith   } else {
282ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
28380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
284ae7cfcebSSatish Balay #else
28571044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
2864ce68768SBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve");
287ae7cfcebSSatish Balay #endif
288da3a660dSBarry Smith   }
2891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2901ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
291b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
293da3a660dSBarry Smith }
2946ee01492SSatish Balay 
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
297dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
298da3a660dSBarry Smith {
299c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
300dfbe8321SBarry Smith   PetscErrorCode ierr;
3014ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
30287828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
303da3a660dSBarry Smith   Vec            tmp = 0;
30467e560aaSBarry Smith 
3053a40ed3dSBarry Smith   PetscFunctionBegin;
3061ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3071ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
308273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
309da3a660dSBarry Smith   if (yy == zz) {
31078b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
311*52e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
31278b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
313da3a660dSBarry Smith   }
31487828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
315752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
316da3a660dSBarry Smith   if (mat->pivots) {
317ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
31880444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
319ae7cfcebSSatish Balay #else
32071044d3cSBarry Smith     LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3212ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
322ae7cfcebSSatish Balay #endif
323a8c6a408SBarry Smith   } else {
324ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
32580444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
326ae7cfcebSSatish Balay #else
32771044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3282ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
329ae7cfcebSSatish Balay #endif
330da3a660dSBarry Smith   }
331a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
332a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
3331ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3341ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
335b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3363a40ed3dSBarry Smith   PetscFunctionReturn(0);
337da3a660dSBarry Smith }
33867e560aaSBarry Smith 
3394a2ae208SSatish Balay #undef __FUNCT__
3404a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
341dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
342da3a660dSBarry Smith {
343c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
3446849ba73SBarry Smith   PetscErrorCode ierr;
3454ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m,one = 1,info;
34687828ca2SBarry Smith   PetscScalar    *x,*y,sone = 1.0;
347da3a660dSBarry Smith   Vec            tmp;
34867e560aaSBarry Smith 
3493a40ed3dSBarry Smith   PetscFunctionBegin;
350273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3511ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3521ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
353da3a660dSBarry Smith   if (yy == zz) {
35478b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
355*52e6d16bSBarry Smith     ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr);
35678b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
357da3a660dSBarry Smith   }
35887828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
359752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
360da3a660dSBarry Smith   if (mat->pivots) {
361ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
36280444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
363ae7cfcebSSatish Balay #else
36471044d3cSBarry Smith     LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info);
3652ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
366ae7cfcebSSatish Balay #endif
3673a40ed3dSBarry Smith   } else {
368ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
36980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
370ae7cfcebSSatish Balay #else
37171044d3cSBarry Smith     LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info);
3722ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
373ae7cfcebSSatish Balay #endif
374da3a660dSBarry Smith   }
37590f02eecSBarry Smith   if (tmp) {
37690f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
37790f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3783a40ed3dSBarry Smith   } else {
37990f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
38090f02eecSBarry Smith   }
3811ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3821ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
383b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3843a40ed3dSBarry Smith   PetscFunctionReturn(0);
385da3a660dSBarry Smith }
386289bc588SBarry Smith /* ------------------------------------------------------------------*/
3874a2ae208SSatish Balay #undef __FUNCT__
3884a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
38913f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
390289bc588SBarry Smith {
391c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
39287828ca2SBarry Smith   PetscScalar    *x,*b,*v = mat->v,zero = 0.0,xt;
393dfbe8321SBarry Smith   PetscErrorCode ierr;
39413f74950SBarry Smith   PetscInt       m = A->m,i;
395aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3964ce68768SBarry Smith   PetscBLASInt   bm = (PetscBLASInt)m, o = 1;
397bc1b551cSSatish Balay #endif
398289bc588SBarry Smith 
3993a40ed3dSBarry Smith   PetscFunctionBegin;
400289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
40171044d3cSBarry Smith     /* this is a hack fix, should have another version without the second BLASdot */
402bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
403289bc588SBarry Smith   }
4041ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4051ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
406b965ef7fSBarry Smith   its  = its*lits;
40777431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
408289bc588SBarry Smith   while (its--) {
409289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
410289bc588SBarry Smith       for (i=0; i<m; i++) {
411aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
412f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
413f1747703SBarry Smith            not happy about returning a double complex */
41413f74950SBarry Smith         PetscInt         _i;
41587828ca2SBarry Smith         PetscScalar sum = b[i];
416f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4173f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
418f1747703SBarry Smith         }
419f1747703SBarry Smith         xt = sum;
420f1747703SBarry Smith #else
42171044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
422f1747703SBarry Smith #endif
42355a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
424289bc588SBarry Smith       }
425289bc588SBarry Smith     }
426289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
427289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
428aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
429f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
430f1747703SBarry Smith            not happy about returning a double complex */
43113f74950SBarry Smith         PetscInt         _i;
43287828ca2SBarry Smith         PetscScalar sum = b[i];
433f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4343f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
435f1747703SBarry Smith         }
436f1747703SBarry Smith         xt = sum;
437f1747703SBarry Smith #else
43871044d3cSBarry Smith         xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o);
439f1747703SBarry Smith #endif
44055a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
441289bc588SBarry Smith       }
442289bc588SBarry Smith     }
443289bc588SBarry Smith   }
4441ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4451ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4463a40ed3dSBarry Smith   PetscFunctionReturn(0);
447289bc588SBarry Smith }
448289bc588SBarry Smith 
449289bc588SBarry Smith /* -----------------------------------------------------------------*/
4504a2ae208SSatish Balay #undef __FUNCT__
4514a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
452dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
453289bc588SBarry Smith {
454c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
45587828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
456dfbe8321SBarry Smith   PetscErrorCode ierr;
4574ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n,_One=1;
458ea709b57SSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
4593a40ed3dSBarry Smith 
4603a40ed3dSBarry Smith   PetscFunctionBegin;
461273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4621ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4631ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
46471044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One);
4651ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4661ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
467b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
469289bc588SBarry Smith }
4706ee01492SSatish Balay 
4714a2ae208SSatish Balay #undef __FUNCT__
4724a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
473dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
474289bc588SBarry Smith {
475c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
47687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
477dfbe8321SBarry Smith   PetscErrorCode ierr;
4784ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
4793a40ed3dSBarry Smith 
4803a40ed3dSBarry Smith   PetscFunctionBegin;
481273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4821ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4831ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
48471044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
487b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->m);
4883a40ed3dSBarry Smith   PetscFunctionReturn(0);
489289bc588SBarry Smith }
4906ee01492SSatish Balay 
4914a2ae208SSatish Balay #undef __FUNCT__
4924a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
493dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
494289bc588SBarry Smith {
495c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
49687828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y,_DOne=1.0;
497dfbe8321SBarry Smith   PetscErrorCode ierr;
4984ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
4993a40ed3dSBarry Smith 
5003a40ed3dSBarry Smith   PetscFunctionBegin;
501273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
502600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5031ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5041ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
50571044d3cSBarry Smith   BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5071ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
508b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
510289bc588SBarry Smith }
5116ee01492SSatish Balay 
5124a2ae208SSatish Balay #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
514dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
515289bc588SBarry Smith {
516c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
51787828ca2SBarry Smith   PetscScalar    *v = mat->v,*x,*y;
518dfbe8321SBarry Smith   PetscErrorCode ierr;
5194ce68768SBarry Smith   PetscBLASInt   m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1;
52087828ca2SBarry Smith   PetscScalar    _DOne=1.0;
5213a40ed3dSBarry Smith 
5223a40ed3dSBarry Smith   PetscFunctionBegin;
523273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
524600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5251ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5261ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
52771044d3cSBarry Smith   BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5281ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5291ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
530b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
5313a40ed3dSBarry Smith   PetscFunctionReturn(0);
532289bc588SBarry Smith }
533289bc588SBarry Smith 
534289bc588SBarry Smith /* -----------------------------------------------------------------*/
5354a2ae208SSatish Balay #undef __FUNCT__
5364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
53713f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
538289bc588SBarry Smith {
539c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
54087828ca2SBarry Smith   PetscScalar    *v;
5416849ba73SBarry Smith   PetscErrorCode ierr;
54213f74950SBarry Smith   PetscInt       i;
54367e560aaSBarry Smith 
5443a40ed3dSBarry Smith   PetscFunctionBegin;
545273d9f13SBarry Smith   *ncols = A->n;
546289bc588SBarry Smith   if (cols) {
54713f74950SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr);
548273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
549289bc588SBarry Smith   }
550289bc588SBarry Smith   if (vals) {
55187828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
552289bc588SBarry Smith     v    = mat->v + row;
5531b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
554289bc588SBarry Smith   }
5553a40ed3dSBarry Smith   PetscFunctionReturn(0);
556289bc588SBarry Smith }
5576ee01492SSatish Balay 
5584a2ae208SSatish Balay #undef __FUNCT__
5594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
56013f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
561289bc588SBarry Smith {
562dfbe8321SBarry Smith   PetscErrorCode ierr;
563606d414cSSatish Balay   PetscFunctionBegin;
564606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
565606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5663a40ed3dSBarry Smith   PetscFunctionReturn(0);
567289bc588SBarry Smith }
568289bc588SBarry Smith /* ----------------------------------------------------------------*/
5694a2ae208SSatish Balay #undef __FUNCT__
5704a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
57113f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
572289bc588SBarry Smith {
573c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
57413f74950SBarry Smith   PetscInt     i,j,idx=0;
575d6dfbf8fSBarry Smith 
5763a40ed3dSBarry Smith   PetscFunctionBegin;
577289bc588SBarry Smith   if (!mat->roworiented) {
578dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
579289bc588SBarry Smith       for (j=0; j<n; j++) {
580cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5812515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58277431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
58358804f6dSBarry Smith #endif
584289bc588SBarry Smith         for (i=0; i<m; i++) {
585cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
5862515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58777431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
58858804f6dSBarry Smith #endif
589cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
590289bc588SBarry Smith         }
591289bc588SBarry Smith       }
5923a40ed3dSBarry Smith     } else {
593289bc588SBarry Smith       for (j=0; j<n; j++) {
594cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
5952515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59677431f27SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
59758804f6dSBarry Smith #endif
598289bc588SBarry Smith         for (i=0; i<m; i++) {
599cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
6002515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60177431f27SBarry Smith           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
60258804f6dSBarry Smith #endif
603cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
604289bc588SBarry Smith         }
605289bc588SBarry Smith       }
606289bc588SBarry Smith     }
6073a40ed3dSBarry Smith   } else {
608dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
609e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
610cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6112515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
61277431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
61358804f6dSBarry Smith #endif
614e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
615cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6162515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
61777431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
61858804f6dSBarry Smith #endif
619cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
620e8d4e0b9SBarry Smith         }
621e8d4e0b9SBarry Smith       }
6223a40ed3dSBarry Smith     } else {
623289bc588SBarry Smith       for (i=0; i<m; i++) {
624cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
6252515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
62677431f27SBarry Smith         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1);
62758804f6dSBarry Smith #endif
628289bc588SBarry Smith         for (j=0; j<n; j++) {
629cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
63177431f27SBarry Smith           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1);
63258804f6dSBarry Smith #endif
633cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
634289bc588SBarry Smith         }
635289bc588SBarry Smith       }
636289bc588SBarry Smith     }
637e8d4e0b9SBarry Smith   }
6383a40ed3dSBarry Smith   PetscFunctionReturn(0);
639289bc588SBarry Smith }
640e8d4e0b9SBarry Smith 
6414a2ae208SSatish Balay #undef __FUNCT__
6424a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
64313f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
644ae80bb75SLois Curfman McInnes {
645ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
64613f74950SBarry Smith   PetscInt     i,j;
64787828ca2SBarry Smith   PetscScalar  *vpt = v;
648ae80bb75SLois Curfman McInnes 
6493a40ed3dSBarry Smith   PetscFunctionBegin;
650ae80bb75SLois Curfman McInnes   /* row-oriented output */
651ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
652ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6531b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
654ae80bb75SLois Curfman McInnes     }
655ae80bb75SLois Curfman McInnes   }
6563a40ed3dSBarry Smith   PetscFunctionReturn(0);
657ae80bb75SLois Curfman McInnes }
658ae80bb75SLois Curfman McInnes 
659289bc588SBarry Smith /* -----------------------------------------------------------------*/
660289bc588SBarry Smith 
661e090d566SSatish Balay #include "petscsys.h"
66288e32edaSLois Curfman McInnes 
6634a2ae208SSatish Balay #undef __FUNCT__
6644a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
665dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A)
66688e32edaSLois Curfman McInnes {
66788e32edaSLois Curfman McInnes   Mat_SeqDense   *a;
66888e32edaSLois Curfman McInnes   Mat            B;
6696849ba73SBarry Smith   PetscErrorCode ierr;
67013f74950SBarry Smith   PetscInt       *scols,i,j,nz,header[4];
67113f74950SBarry Smith   int            fd;
67213f74950SBarry Smith   PetscMPIInt    size;
67313f74950SBarry Smith   PetscInt       *rowlengths = 0,M,N,*cols;
67487828ca2SBarry Smith   PetscScalar    *vals,*svals,*v,*w;
67519bcc07fSBarry Smith   MPI_Comm       comm = ((PetscObject)viewer)->comm;
67688e32edaSLois Curfman McInnes 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
678d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
67929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
680b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6810752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
682552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
68388e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
68488e32edaSLois Curfman McInnes 
68590ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
6865c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
687be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6885c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
68990ace30eSBarry Smith     B    = *A;
69090ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6914a41a4d3SSatish Balay     v    = a->v;
6924a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6934a41a4d3SSatish Balay        from row major to column major */
694b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
69590ace30eSBarry Smith     /* read in nonzero values */
6964a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6974a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6984a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6994a41a4d3SSatish Balay       for (i=0; i<M; i++) {
7004a41a4d3SSatish Balay         *v++ =w[i*N+j];
7014a41a4d3SSatish Balay       }
7024a41a4d3SSatish Balay     }
703606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
7046d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7056d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70690ace30eSBarry Smith   } else {
70788e32edaSLois Curfman McInnes     /* read row lengths */
70813f74950SBarry Smith     ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
7090752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
71088e32edaSLois Curfman McInnes 
71188e32edaSLois Curfman McInnes     /* create our matrix */
7125c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
713be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
7145c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
71588e32edaSLois Curfman McInnes     B = *A;
71688e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
71788e32edaSLois Curfman McInnes     v = a->v;
71888e32edaSLois Curfman McInnes 
71988e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
72013f74950SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr);
721b0a32e0cSBarry Smith     cols = scols;
7220752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
72387828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
724b0a32e0cSBarry Smith     vals = svals;
7250752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
72688e32edaSLois Curfman McInnes 
72788e32edaSLois Curfman McInnes     /* insert into matrix */
72888e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
72988e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
73088e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
73188e32edaSLois Curfman McInnes     }
732606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
733606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
734606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
73588e32edaSLois Curfman McInnes 
7366d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7376d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73890ace30eSBarry Smith   }
7393a40ed3dSBarry Smith   PetscFunctionReturn(0);
74088e32edaSLois Curfman McInnes }
74188e32edaSLois Curfman McInnes 
742e090d566SSatish Balay #include "petscsys.h"
743932b0c3eSLois Curfman McInnes 
7444a2ae208SSatish Balay #undef __FUNCT__
7454a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
7466849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
747289bc588SBarry Smith {
748932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
749dfbe8321SBarry Smith   PetscErrorCode    ierr;
75013f74950SBarry Smith   PetscInt          i,j;
751fb9695e5SSatish Balay   char              *name;
75287828ca2SBarry Smith   PetscScalar       *v;
753f3ef73ceSBarry Smith   PetscViewerFormat format;
754932b0c3eSLois Curfman McInnes 
7553a40ed3dSBarry Smith   PetscFunctionBegin;
756435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
757b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
758456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7593a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
760fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
761b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
762273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
76344cd7ae7SLois Curfman McInnes       v = a->v + i;
76477431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
765273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
766aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
767329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
76877431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
769329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
77077431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7716831982aSBarry Smith         }
77280cd9d93SLois Curfman McInnes #else
7736831982aSBarry Smith         if (*v) {
77477431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,*v);CHKERRQ(ierr);
7756831982aSBarry Smith         }
77680cd9d93SLois Curfman McInnes #endif
7771b807ce4Svictorle         v += a->lda;
77880cd9d93SLois Curfman McInnes       }
779b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
78080cd9d93SLois Curfman McInnes     }
781b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7823a40ed3dSBarry Smith   } else {
783b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
784aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
785ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
78647989497SBarry Smith     /* determine if matrix has all real values */
78747989497SBarry Smith     v = a->v;
788201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
789ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
79047989497SBarry Smith     }
79147989497SBarry Smith #endif
792fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7933a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
79477431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->m,A->n);CHKERRQ(ierr);
79577431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->m,A->n);CHKERRQ(ierr);
796fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
797ffac6cdbSBarry Smith     }
798ffac6cdbSBarry Smith 
799273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
800932b0c3eSLois Curfman McInnes       v = a->v + i;
801273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
802aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
80347989497SBarry Smith         if (allreal) {
8041b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
80547989497SBarry Smith         } else {
8061b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
80747989497SBarry Smith         }
808289bc588SBarry Smith #else
8091b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
810289bc588SBarry Smith #endif
8111b807ce4Svictorle 	v += a->lda;
812289bc588SBarry Smith       }
813b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
814289bc588SBarry Smith     }
815fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
816b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
817ffac6cdbSBarry Smith     }
818b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
819da3a660dSBarry Smith   }
820b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8213a40ed3dSBarry Smith   PetscFunctionReturn(0);
822289bc588SBarry Smith }
823289bc588SBarry Smith 
8244a2ae208SSatish Balay #undef __FUNCT__
8254a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
8266849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
827932b0c3eSLois Curfman McInnes {
828932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
8296849ba73SBarry Smith   PetscErrorCode    ierr;
83013f74950SBarry Smith   int               fd;
83113f74950SBarry Smith   PetscInt          ict,j,n = A->n,m = A->m,i,*col_lens,nz = m*n;
83287828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
833f3ef73ceSBarry Smith   PetscViewerFormat format;
834932b0c3eSLois Curfman McInnes 
8353a40ed3dSBarry Smith   PetscFunctionBegin;
836b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
83790ace30eSBarry Smith 
838b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
839fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
84090ace30eSBarry Smith     /* store the matrix as a dense matrix */
84113f74950SBarry Smith     ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
842552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
84390ace30eSBarry Smith     col_lens[1] = m;
84490ace30eSBarry Smith     col_lens[2] = n;
84590ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8466f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
847606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
84890ace30eSBarry Smith 
84990ace30eSBarry Smith     /* write out matrix, by rows */
85087828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
85190ace30eSBarry Smith     v    = a->v;
85290ace30eSBarry Smith     for (i=0; i<m; i++) {
85390ace30eSBarry Smith       for (j=0; j<n; j++) {
85490ace30eSBarry Smith         vals[i + j*m] = *v++;
85590ace30eSBarry Smith       }
85690ace30eSBarry Smith     }
8576f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
858606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
85990ace30eSBarry Smith   } else {
86013f74950SBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr);
861552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
862932b0c3eSLois Curfman McInnes     col_lens[1] = m;
863932b0c3eSLois Curfman McInnes     col_lens[2] = n;
864932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
865932b0c3eSLois Curfman McInnes 
866932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
867932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8686f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
869932b0c3eSLois Curfman McInnes 
870932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
871932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
872932b0c3eSLois Curfman McInnes     ict = 0;
873932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
874932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
875932b0c3eSLois Curfman McInnes     }
8766f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
877606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
878932b0c3eSLois Curfman McInnes 
879932b0c3eSLois Curfman McInnes     /* store nonzero values */
88087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
881932b0c3eSLois Curfman McInnes     ict  = 0;
882932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
883932b0c3eSLois Curfman McInnes       v = a->v + i;
884932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8851b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
886932b0c3eSLois Curfman McInnes       }
887932b0c3eSLois Curfman McInnes     }
8886f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
889606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
89090ace30eSBarry Smith   }
8913a40ed3dSBarry Smith   PetscFunctionReturn(0);
892932b0c3eSLois Curfman McInnes }
893932b0c3eSLois Curfman McInnes 
8944a2ae208SSatish Balay #undef __FUNCT__
8954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
896dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
897f1af5d2fSBarry Smith {
898f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
899f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
9006849ba73SBarry Smith   PetscErrorCode    ierr;
90113f74950SBarry Smith   PetscInt          m = A->m,n = A->n,color,i,j;
90287828ca2SBarry Smith   PetscScalar       *v = a->v;
903b0a32e0cSBarry Smith   PetscViewer       viewer;
904b0a32e0cSBarry Smith   PetscDraw         popup;
905329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
906f3ef73ceSBarry Smith   PetscViewerFormat format;
907f1af5d2fSBarry Smith 
908f1af5d2fSBarry Smith   PetscFunctionBegin;
909f1af5d2fSBarry Smith 
910f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
911b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
912b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
913f1af5d2fSBarry Smith 
914f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
915fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
916f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
917b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
918f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
919f1af5d2fSBarry Smith       x_l = j;
920f1af5d2fSBarry Smith       x_r = x_l + 1.0;
921f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
922f1af5d2fSBarry Smith         y_l = m - i - 1.0;
923f1af5d2fSBarry Smith         y_r = y_l + 1.0;
924f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
925329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
926b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
927329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
928b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
929f1af5d2fSBarry Smith         } else {
930f1af5d2fSBarry Smith           continue;
931f1af5d2fSBarry Smith         }
932f1af5d2fSBarry Smith #else
933f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
934b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
935f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
936b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
937f1af5d2fSBarry Smith         } else {
938f1af5d2fSBarry Smith           continue;
939f1af5d2fSBarry Smith         }
940f1af5d2fSBarry Smith #endif
941b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
942f1af5d2fSBarry Smith       }
943f1af5d2fSBarry Smith     }
944f1af5d2fSBarry Smith   } else {
945f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
946f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
947f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
948f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
949f1af5d2fSBarry Smith     }
950b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
951b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
952b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
953f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
954f1af5d2fSBarry Smith       x_l = j;
955f1af5d2fSBarry Smith       x_r = x_l + 1.0;
956f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
957f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
958f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
959b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
960b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
961f1af5d2fSBarry Smith       }
962f1af5d2fSBarry Smith     }
963f1af5d2fSBarry Smith   }
964f1af5d2fSBarry Smith   PetscFunctionReturn(0);
965f1af5d2fSBarry Smith }
966f1af5d2fSBarry Smith 
9674a2ae208SSatish Balay #undef __FUNCT__
9684a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
969dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
970f1af5d2fSBarry Smith {
971b0a32e0cSBarry Smith   PetscDraw      draw;
972f1af5d2fSBarry Smith   PetscTruth     isnull;
973329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
974dfbe8321SBarry Smith   PetscErrorCode ierr;
975f1af5d2fSBarry Smith 
976f1af5d2fSBarry Smith   PetscFunctionBegin;
977b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
978b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
979abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
980f1af5d2fSBarry Smith 
981f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
982273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
983f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
984b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
985b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
986f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
987f1af5d2fSBarry Smith   PetscFunctionReturn(0);
988f1af5d2fSBarry Smith }
989f1af5d2fSBarry Smith 
9904a2ae208SSatish Balay #undef __FUNCT__
9914a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
992dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
993932b0c3eSLois Curfman McInnes {
994932b0c3eSLois Curfman McInnes   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
995dfbe8321SBarry Smith   PetscErrorCode ierr;
99632077d6dSBarry Smith   PetscTruth     issocket,iascii,isbinary,isdraw;
997932b0c3eSLois Curfman McInnes 
9983a40ed3dSBarry Smith   PetscFunctionBegin;
999b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
100032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1001fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
1002fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
10030f5bd95cSBarry Smith 
10040f5bd95cSBarry Smith   if (issocket) {
1005634064b4SBarry Smith     if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
1006b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
100732077d6dSBarry Smith   } else if (iascii) {
10083a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
10090f5bd95cSBarry Smith   } else if (isbinary) {
10103a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1011f1af5d2fSBarry Smith   } else if (isdraw) {
1012f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
10135cd90555SBarry Smith   } else {
1014958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
1015932b0c3eSLois Curfman McInnes   }
10163a40ed3dSBarry Smith   PetscFunctionReturn(0);
1017932b0c3eSLois Curfman McInnes }
1018289bc588SBarry Smith 
10194a2ae208SSatish Balay #undef __FUNCT__
10204a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1021dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat)
1022289bc588SBarry Smith {
1023ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1024dfbe8321SBarry Smith   PetscErrorCode ierr;
102590f02eecSBarry Smith 
10263a40ed3dSBarry Smith   PetscFunctionBegin;
1027aa482453SBarry Smith #if defined(PETSC_USE_LOG)
102877431f27SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->m,mat->n);
1029a5a9c739SBarry Smith #endif
1030606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1031606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1032606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
1033901853e0SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
10343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1035289bc588SBarry Smith }
1036289bc588SBarry Smith 
10374a2ae208SSatish Balay #undef __FUNCT__
10384a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
1039dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout)
1040289bc588SBarry Smith {
1041c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
10426849ba73SBarry Smith   PetscErrorCode ierr;
104313f74950SBarry Smith   PetscInt       k,j,m,n,M;
104487828ca2SBarry Smith   PetscScalar    *v,tmp;
104548b35521SBarry Smith 
10463a40ed3dSBarry Smith   PetscFunctionBegin;
10471b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10487c922b88SBarry Smith   if (!matout) { /* in place transpose */
1049a5ce6ee0Svictorle     if (m != n) {
1050634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1051d64ed03dSBarry Smith     } else {
1052d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1053289bc588SBarry Smith         for (k=0; k<j; k++) {
10541b807ce4Svictorle           tmp = v[j + k*M];
10551b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10561b807ce4Svictorle           v[k + j*M] = tmp;
1057289bc588SBarry Smith         }
1058289bc588SBarry Smith       }
1059d64ed03dSBarry Smith     }
10603a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1061d3e5ee88SLois Curfman McInnes     Mat          tmat;
1062ec8511deSBarry Smith     Mat_SeqDense *tmatd;
106387828ca2SBarry Smith     PetscScalar  *v2;
1064ea709b57SSatish Balay 
10655c5985e7SKris Buschelman     ierr  = MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);CHKERRQ(ierr);
10665c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10675c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1068ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10690de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1070d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10711b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1072d3e5ee88SLois Curfman McInnes     }
10736d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10746d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1075d3e5ee88SLois Curfman McInnes     *matout = tmat;
107648b35521SBarry Smith   }
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078289bc588SBarry Smith }
1079289bc588SBarry Smith 
10804a2ae208SSatish Balay #undef __FUNCT__
10814a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
1082dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1083289bc588SBarry Smith {
1084c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1085c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
108613f74950SBarry Smith   PetscInt     i,j;
108787828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10889ea5d5aeSSatish Balay 
10893a40ed3dSBarry Smith   PetscFunctionBegin;
1090273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1091273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10921b807ce4Svictorle   for (i=0; i<A1->m; i++) {
10931b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
10941b807ce4Svictorle     for (j=0; j<A1->n; j++) {
10953a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10961b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
10971b807ce4Svictorle     }
1098289bc588SBarry Smith   }
109977c4ece6SBarry Smith   *flg = PETSC_TRUE;
11003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1101289bc588SBarry Smith }
1102289bc588SBarry Smith 
11034a2ae208SSatish Balay #undef __FUNCT__
11044a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
1105dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1106289bc588SBarry Smith {
1107c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1108dfbe8321SBarry Smith   PetscErrorCode ierr;
110913f74950SBarry Smith   PetscInt       i,n,len;
111087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
111144cd7ae7SLois Curfman McInnes 
11123a40ed3dSBarry Smith   PetscFunctionBegin;
11137a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
11147a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
11151ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1116273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1117273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
111844cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
11191b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1120289bc588SBarry Smith   }
11211ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
11223a40ed3dSBarry Smith   PetscFunctionReturn(0);
1123289bc588SBarry Smith }
1124289bc588SBarry Smith 
11254a2ae208SSatish Balay #undef __FUNCT__
11264a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
1127dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1128289bc588SBarry Smith {
1129c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
113087828ca2SBarry Smith   PetscScalar    *l,*r,x,*v;
1131dfbe8321SBarry Smith   PetscErrorCode ierr;
113213f74950SBarry Smith   PetscInt       i,j,m = A->m,n = A->n;
113355659b69SBarry Smith 
11343a40ed3dSBarry Smith   PetscFunctionBegin;
113528988994SBarry Smith   if (ll) {
11367a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11371ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1138273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1139da3a660dSBarry Smith     for (i=0; i<m; i++) {
1140da3a660dSBarry Smith       x = l[i];
1141da3a660dSBarry Smith       v = mat->v + i;
1142da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1143da3a660dSBarry Smith     }
11441ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1145b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1146da3a660dSBarry Smith   }
114728988994SBarry Smith   if (rr) {
11487a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11491ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1150273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1151da3a660dSBarry Smith     for (i=0; i<n; i++) {
1152da3a660dSBarry Smith       x = r[i];
1153da3a660dSBarry Smith       v = mat->v + i*m;
1154da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1155da3a660dSBarry Smith     }
11561ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1157b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1158da3a660dSBarry Smith   }
11593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1160289bc588SBarry Smith }
1161289bc588SBarry Smith 
11624a2ae208SSatish Balay #undef __FUNCT__
11634a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1164dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1165289bc588SBarry Smith {
1166c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
116787828ca2SBarry Smith   PetscScalar  *v = mat->v;
1168329f5518SBarry Smith   PetscReal    sum = 0.0;
116913f74950SBarry Smith   PetscInt     lda=mat->lda,m=A->m,i,j;
117055659b69SBarry Smith 
11713a40ed3dSBarry Smith   PetscFunctionBegin;
1172289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1173a5ce6ee0Svictorle     if (lda>m) {
1174a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1175a5ce6ee0Svictorle 	v = mat->v+j*lda;
1176a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1177a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1178a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1179a5ce6ee0Svictorle #else
1180a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1181a5ce6ee0Svictorle #endif
1182a5ce6ee0Svictorle 	}
1183a5ce6ee0Svictorle       }
1184a5ce6ee0Svictorle     } else {
1185273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1186aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1187329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1188289bc588SBarry Smith #else
1189289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1190289bc588SBarry Smith #endif
1191289bc588SBarry Smith       }
1192a5ce6ee0Svictorle     }
1193064f8208SBarry Smith     *nrm = sqrt(sum);
1194b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
11953a40ed3dSBarry Smith   } else if (type == NORM_1) {
1196064f8208SBarry Smith     *nrm = 0.0;
1197273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
11981b807ce4Svictorle       v = mat->v + j*mat->lda;
1199289bc588SBarry Smith       sum = 0.0;
1200273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
120133a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1202289bc588SBarry Smith       }
1203064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1204289bc588SBarry Smith     }
1205b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
12063a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1207064f8208SBarry Smith     *nrm = 0.0;
1208273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1209289bc588SBarry Smith       v = mat->v + j;
1210289bc588SBarry Smith       sum = 0.0;
1211273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
12121b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1213289bc588SBarry Smith       }
1214064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1215289bc588SBarry Smith     }
1216b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
12173a40ed3dSBarry Smith   } else {
121829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1219289bc588SBarry Smith   }
12203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1221289bc588SBarry Smith }
1222289bc588SBarry Smith 
12234a2ae208SSatish Balay #undef __FUNCT__
12244a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
1225dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op)
1226289bc588SBarry Smith {
1227c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
122867e560aaSBarry Smith 
12293a40ed3dSBarry Smith   PetscFunctionBegin;
1230b5a2b587SKris Buschelman   switch (op) {
1231b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1232b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1233b5a2b587SKris Buschelman     break;
1234b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1235b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1236b5a2b587SKris Buschelman     break;
1237b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1238b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1239b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1240b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1241b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1242b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1243b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1244b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1245b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1246b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1247b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1248b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1249b5a2b587SKris Buschelman     break;
125077e54ba9SKris Buschelman   case MAT_SYMMETRIC:
125177e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12529a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12539a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12549a4540c5SBarry Smith   case MAT_HERMITIAN:
12559a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12569a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12579a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
125877e54ba9SKris Buschelman     break;
1259b5a2b587SKris Buschelman   default:
126029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12613a40ed3dSBarry Smith   }
12623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1263289bc588SBarry Smith }
1264289bc588SBarry Smith 
12654a2ae208SSatish Balay #undef __FUNCT__
12664a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
1267dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A)
12686f0a148fSBarry Smith {
1269ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12706849ba73SBarry Smith   PetscErrorCode ierr;
127113f74950SBarry Smith   PetscInt       lda=l->lda,m=A->m,j;
12723a40ed3dSBarry Smith 
12733a40ed3dSBarry Smith   PetscFunctionBegin;
1274a5ce6ee0Svictorle   if (lda>m) {
1275a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1276a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1277a5ce6ee0Svictorle     }
1278a5ce6ee0Svictorle   } else {
127987828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1280a5ce6ee0Svictorle   }
12813a40ed3dSBarry Smith   PetscFunctionReturn(0);
12826f0a148fSBarry Smith }
12836f0a148fSBarry Smith 
12844a2ae208SSatish Balay #undef __FUNCT__
12854a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1286dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
12876f0a148fSBarry Smith {
1288ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
12896849ba73SBarry Smith   PetscErrorCode ierr;
129013f74950SBarry Smith   PetscInt       n = A->n,i,j,N,*rows;
129187828ca2SBarry Smith   PetscScalar    *slot;
129255659b69SBarry Smith 
12933a40ed3dSBarry Smith   PetscFunctionBegin;
1294e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
129578b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
12966f0a148fSBarry Smith   for (i=0; i<N; i++) {
12976f0a148fSBarry Smith     slot = l->v + rows[i];
12986f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12996f0a148fSBarry Smith   }
13006f0a148fSBarry Smith   if (diag) {
13016f0a148fSBarry Smith     for (i=0; i<N; i++) {
13026f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1303260925b8SBarry Smith       *slot = *diag;
13046f0a148fSBarry Smith     }
13056f0a148fSBarry Smith   }
1306260925b8SBarry Smith   ISRestoreIndices(is,&rows);
13073a40ed3dSBarry Smith   PetscFunctionReturn(0);
13086f0a148fSBarry Smith }
1309557bce09SLois Curfman McInnes 
13104a2ae208SSatish Balay #undef __FUNCT__
13114a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
1312dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[])
131364e87e97SBarry Smith {
1314c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
13153a40ed3dSBarry Smith 
13163a40ed3dSBarry Smith   PetscFunctionBegin;
131764e87e97SBarry Smith   *array = mat->v;
13183a40ed3dSBarry Smith   PetscFunctionReturn(0);
131964e87e97SBarry Smith }
13200754003eSLois Curfman McInnes 
13214a2ae208SSatish Balay #undef __FUNCT__
13224a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
1323dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1324ff14e315SSatish Balay {
13253a40ed3dSBarry Smith   PetscFunctionBegin;
132609b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1328ff14e315SSatish Balay }
13290754003eSLois Curfman McInnes 
13304a2ae208SSatish Balay #undef __FUNCT__
13314a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
133213f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
13330754003eSLois Curfman McInnes {
1334c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
13356849ba73SBarry Smith   PetscErrorCode ierr;
133613f74950SBarry Smith   PetscInt       i,j,m = A->m,*irow,*icol,nrows,ncols;
133787828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
13380754003eSLois Curfman McInnes   Mat            newmat;
13390754003eSLois Curfman McInnes 
13403a40ed3dSBarry Smith   PetscFunctionBegin;
134178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
134278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1343e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1344e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13450754003eSLois Curfman McInnes 
1346182d2002SSatish Balay   /* Check submatrixcall */
1347182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
134813f74950SBarry Smith     PetscInt n_cols,n_rows;
1349182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
135029bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1351182d2002SSatish Balay     newmat = *B;
1352182d2002SSatish Balay   } else {
13530754003eSLois Curfman McInnes     /* Create and fill new matrix */
13545c5985e7SKris Buschelman     ierr = MatCreate(A->comm,nrows,ncols,nrows,ncols,&newmat);CHKERRQ(ierr);
13555c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13565c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1357182d2002SSatish Balay   }
1358182d2002SSatish Balay 
1359182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1360182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1361182d2002SSatish Balay 
1362182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1363182d2002SSatish Balay     av = v + m*icol[i];
1364182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1365182d2002SSatish Balay       *bv++ = av[irow[j]];
13660754003eSLois Curfman McInnes     }
13670754003eSLois Curfman McInnes   }
1368182d2002SSatish Balay 
1369182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13706d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13716d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13720754003eSLois Curfman McInnes 
13730754003eSLois Curfman McInnes   /* Free work space */
137478b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
137578b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1376182d2002SSatish Balay   *B = newmat;
13773a40ed3dSBarry Smith   PetscFunctionReturn(0);
13780754003eSLois Curfman McInnes }
13790754003eSLois Curfman McInnes 
13804a2ae208SSatish Balay #undef __FUNCT__
13814a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
138213f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1383905e6a2fSBarry Smith {
13846849ba73SBarry Smith   PetscErrorCode ierr;
138513f74950SBarry Smith   PetscInt       i;
1386905e6a2fSBarry Smith 
13873a40ed3dSBarry Smith   PetscFunctionBegin;
1388905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1389b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1390905e6a2fSBarry Smith   }
1391905e6a2fSBarry Smith 
1392905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13936a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1394905e6a2fSBarry Smith   }
13953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1396905e6a2fSBarry Smith }
1397905e6a2fSBarry Smith 
13984a2ae208SSatish Balay #undef __FUNCT__
13994a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1400dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
14014b0e389bSBarry Smith {
14024b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
14036849ba73SBarry Smith   PetscErrorCode ierr;
140413f74950SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j;
14053a40ed3dSBarry Smith 
14063a40ed3dSBarry Smith   PetscFunctionBegin;
140733f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
140833f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1409cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
14103a40ed3dSBarry Smith     PetscFunctionReturn(0);
14113a40ed3dSBarry Smith   }
14120dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1413a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
14140dbb7854Svictorle     for (j=0; j<n; j++) {
1415a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1416a5ce6ee0Svictorle     }
1417a5ce6ee0Svictorle   } else {
141887828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1419a5ce6ee0Svictorle   }
1420273d9f13SBarry Smith   PetscFunctionReturn(0);
1421273d9f13SBarry Smith }
1422273d9f13SBarry Smith 
14234a2ae208SSatish Balay #undef __FUNCT__
14244a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1425dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A)
1426273d9f13SBarry Smith {
1427dfbe8321SBarry Smith   PetscErrorCode ierr;
1428273d9f13SBarry Smith 
1429273d9f13SBarry Smith   PetscFunctionBegin;
1430273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14313a40ed3dSBarry Smith   PetscFunctionReturn(0);
14324b0e389bSBarry Smith }
14334b0e389bSBarry Smith 
1434289bc588SBarry Smith /* -------------------------------------------------------------------*/
1435a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1436905e6a2fSBarry Smith        MatGetRow_SeqDense,
1437905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1438905e6a2fSBarry Smith        MatMult_SeqDense,
143997304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14407c922b88SBarry Smith        MatMultTranspose_SeqDense,
14417c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1442905e6a2fSBarry Smith        MatSolve_SeqDense,
1443905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14447c922b88SBarry Smith        MatSolveTranspose_SeqDense,
144597304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1446905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1447905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1448ec8511deSBarry Smith        MatRelax_SeqDense,
1449ec8511deSBarry Smith        MatTranspose_SeqDense,
145097304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1451905e6a2fSBarry Smith        MatEqual_SeqDense,
1452905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1453905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1454905e6a2fSBarry Smith        MatNorm_SeqDense,
145597304618SKris Buschelman /*20*/ 0,
1456905e6a2fSBarry Smith        0,
1457905e6a2fSBarry Smith        0,
1458905e6a2fSBarry Smith        MatSetOption_SeqDense,
1459905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
146097304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1461905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1462905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1463905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1464905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
146597304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1466273d9f13SBarry Smith        0,
1467905e6a2fSBarry Smith        0,
1468905e6a2fSBarry Smith        MatGetArray_SeqDense,
1469905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
147097304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1471a5ae1ecdSBarry Smith        0,
1472a5ae1ecdSBarry Smith        0,
1473a5ae1ecdSBarry Smith        0,
1474a5ae1ecdSBarry Smith        0,
147597304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1476a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1477a5ae1ecdSBarry Smith        0,
14784b0e389bSBarry Smith        MatGetValues_SeqDense,
1479a5ae1ecdSBarry Smith        MatCopy_SeqDense,
148097304618SKris Buschelman /*45*/ 0,
1481a5ae1ecdSBarry Smith        MatScale_SeqDense,
1482a5ae1ecdSBarry Smith        0,
1483a5ae1ecdSBarry Smith        0,
1484a5ae1ecdSBarry Smith        0,
1485521d7252SBarry Smith /*50*/ 0,
1486a5ae1ecdSBarry Smith        0,
1487a5ae1ecdSBarry Smith        0,
1488a5ae1ecdSBarry Smith        0,
1489a5ae1ecdSBarry Smith        0,
149097304618SKris Buschelman /*55*/ 0,
1491a5ae1ecdSBarry Smith        0,
1492a5ae1ecdSBarry Smith        0,
1493a5ae1ecdSBarry Smith        0,
1494a5ae1ecdSBarry Smith        0,
149597304618SKris Buschelman /*60*/ 0,
1496e03a110bSBarry Smith        MatDestroy_SeqDense,
1497e03a110bSBarry Smith        MatView_SeqDense,
149897304618SKris Buschelman        MatGetPetscMaps_Petsc,
149997304618SKris Buschelman        0,
150097304618SKris Buschelman /*65*/ 0,
150197304618SKris Buschelman        0,
150297304618SKris Buschelman        0,
150397304618SKris Buschelman        0,
150497304618SKris Buschelman        0,
150597304618SKris Buschelman /*70*/ 0,
150697304618SKris Buschelman        0,
150797304618SKris Buschelman        0,
150897304618SKris Buschelman        0,
150997304618SKris Buschelman        0,
151097304618SKris Buschelman /*75*/ 0,
151197304618SKris Buschelman        0,
151297304618SKris Buschelman        0,
151397304618SKris Buschelman        0,
151497304618SKris Buschelman        0,
151597304618SKris Buschelman /*80*/ 0,
151697304618SKris Buschelman        0,
151797304618SKris Buschelman        0,
151897304618SKris Buschelman        0,
1519865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense,
1520865e5f61SKris Buschelman        0,
1521865e5f61SKris Buschelman        0,
1522865e5f61SKris Buschelman        0,
1523865e5f61SKris Buschelman        0,
1524865e5f61SKris Buschelman        0,
1525865e5f61SKris Buschelman /*90*/ 0,
1526865e5f61SKris Buschelman        0,
1527865e5f61SKris Buschelman        0,
1528865e5f61SKris Buschelman        0,
1529865e5f61SKris Buschelman        0,
1530865e5f61SKris Buschelman /*95*/ 0,
1531865e5f61SKris Buschelman        0,
1532865e5f61SKris Buschelman        0,
1533865e5f61SKris Buschelman        0};
153490ace30eSBarry Smith 
15354a2ae208SSatish Balay #undef __FUNCT__
15364a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
15374b828684SBarry Smith /*@C
1538fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1539d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1540d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1541289bc588SBarry Smith 
1542db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1543db81eaa0SLois Curfman McInnes 
154420563c6bSBarry Smith    Input Parameters:
1545db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15460c775827SLois Curfman McInnes .  m - number of rows
154718f449edSLois Curfman McInnes .  n - number of columns
1548db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1549dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
155020563c6bSBarry Smith 
155120563c6bSBarry Smith    Output Parameter:
155244cd7ae7SLois Curfman McInnes .  A - the matrix
155320563c6bSBarry Smith 
1554b259b22eSLois Curfman McInnes    Notes:
155518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
155618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1557b4fd4287SBarry Smith    set data=PETSC_NULL.
155818f449edSLois Curfman McInnes 
1559027ccd11SLois Curfman McInnes    Level: intermediate
1560027ccd11SLois Curfman McInnes 
1561dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1562d65003e9SLois Curfman McInnes 
1563db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
156420563c6bSBarry Smith @*/
156513f74950SBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
1566289bc588SBarry Smith {
1567dfbe8321SBarry Smith   PetscErrorCode ierr;
15683b2fbd54SBarry Smith 
15693a40ed3dSBarry Smith   PetscFunctionBegin;
1570273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1571273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1572273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1573273d9f13SBarry Smith   PetscFunctionReturn(0);
1574273d9f13SBarry Smith }
1575273d9f13SBarry Smith 
15764a2ae208SSatish Balay #undef __FUNCT__
15774a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1578273d9f13SBarry Smith /*@C
1579273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1580273d9f13SBarry Smith 
1581273d9f13SBarry Smith    Collective on MPI_Comm
1582273d9f13SBarry Smith 
1583273d9f13SBarry Smith    Input Parameters:
1584273d9f13SBarry Smith +  A - the matrix
1585273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1586273d9f13SBarry Smith 
1587273d9f13SBarry Smith    Notes:
1588273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1589273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1590273d9f13SBarry Smith    set data=PETSC_NULL.
1591273d9f13SBarry Smith 
1592273d9f13SBarry Smith    Level: intermediate
1593273d9f13SBarry Smith 
1594273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1595273d9f13SBarry Smith 
1596273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1597273d9f13SBarry Smith @*/
1598dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1599273d9f13SBarry Smith {
1600dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscScalar[]);
1601a23d5eceSKris Buschelman 
1602a23d5eceSKris Buschelman   PetscFunctionBegin;
1603a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1604a23d5eceSKris Buschelman   if (f) {
1605a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1606a23d5eceSKris Buschelman   }
1607a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1608a23d5eceSKris Buschelman }
1609a23d5eceSKris Buschelman 
1610a23d5eceSKris Buschelman EXTERN_C_BEGIN
1611a23d5eceSKris Buschelman #undef __FUNCT__
1612a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1613dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1614a23d5eceSKris Buschelman {
1615273d9f13SBarry Smith   Mat_SeqDense   *b;
1616dfbe8321SBarry Smith   PetscErrorCode ierr;
1617273d9f13SBarry Smith 
1618273d9f13SBarry Smith   PetscFunctionBegin;
1619273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1620273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1621273d9f13SBarry Smith   if (!data) {
162287828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
162387828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1624273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
1625*52e6d16bSBarry Smith     ierr = PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));CHKERRQ(ierr);
1626273d9f13SBarry Smith   } else { /* user-allocated storage */
1627273d9f13SBarry Smith     b->v          = data;
1628273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1629273d9f13SBarry Smith   }
1630273d9f13SBarry Smith   PetscFunctionReturn(0);
1631273d9f13SBarry Smith }
1632a23d5eceSKris Buschelman EXTERN_C_END
1633273d9f13SBarry Smith 
16341b807ce4Svictorle #undef __FUNCT__
16351b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
16361b807ce4Svictorle /*@C
16371b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
16381b807ce4Svictorle 
16391b807ce4Svictorle   Input parameter:
16401b807ce4Svictorle + A - the matrix
16411b807ce4Svictorle - lda - the leading dimension
16421b807ce4Svictorle 
16431b807ce4Svictorle   Notes:
16441b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16451b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16461b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16471b807ce4Svictorle 
16481b807ce4Svictorle   Level: intermediate
16491b807ce4Svictorle 
16501b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16511b807ce4Svictorle 
16521b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16531b807ce4Svictorle @*/
165413f74950SBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda)
16551b807ce4Svictorle {
16561b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16571b807ce4Svictorle   PetscFunctionBegin;
165877431f27SBarry Smith   if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->m);
16591b807ce4Svictorle   b->lda = lda;
16601b807ce4Svictorle   PetscFunctionReturn(0);
16611b807ce4Svictorle }
16621b807ce4Svictorle 
16630bad9183SKris Buschelman /*MC
1664fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16650bad9183SKris Buschelman 
16660bad9183SKris Buschelman    Options Database Keys:
16670bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16680bad9183SKris Buschelman 
16690bad9183SKris Buschelman   Level: beginner
16700bad9183SKris Buschelman 
16710bad9183SKris Buschelman .seealso: MatCreateSeqDense
16720bad9183SKris Buschelman M*/
16730bad9183SKris Buschelman 
1674273d9f13SBarry Smith EXTERN_C_BEGIN
16754a2ae208SSatish Balay #undef __FUNCT__
16764a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1677dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B)
1678273d9f13SBarry Smith {
1679273d9f13SBarry Smith   Mat_SeqDense   *b;
1680dfbe8321SBarry Smith   PetscErrorCode ierr;
16817c334f02SBarry Smith   PetscMPIInt    size;
1682273d9f13SBarry Smith 
1683273d9f13SBarry Smith   PetscFunctionBegin;
1684273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
168529bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
168655659b69SBarry Smith 
1687273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1688273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1689273d9f13SBarry Smith 
1690b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1691549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
169244cd7ae7SLois Curfman McInnes   B->factor       = 0;
169390f02eecSBarry Smith   B->mapping      = 0;
1694*52e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr);
169544cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
169618f449edSLois Curfman McInnes 
16978a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
16988a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1699a5ae1ecdSBarry Smith 
170044cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1701273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1702273d9f13SBarry Smith   b->v            = 0;
17031b807ce4Svictorle   b->lda          = B->m;
17044e220ebcSLois Curfman McInnes 
1705a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1706a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1707a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
17083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1709289bc588SBarry Smith }
1710273d9f13SBarry Smith EXTERN_C_END
1711