xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 634064b45b5c838063ae82f97ffb7e99245dcdb5)
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"
10268466fbSBarry Smith int 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;
13a5ce6ee0Svictorle   int          N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1;
143a40ed3dSBarry Smith 
153a40ed3dSBarry Smith   PetscFunctionBegin;
16a5ce6ee0Svictorle   if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
17a5ce6ee0Svictorle   if (ldax>m || lday>m) {
18a5ce6ee0Svictorle     for (j=0; j<X->n; j++) {
19268466fbSBarry Smith       BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one);
20a5ce6ee0Svictorle     }
21a5ce6ee0Svictorle   } else {
22268466fbSBarry Smith     BLaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one);
23a5ce6ee0Svictorle   }
24b0a32e0cSBarry Smith   PetscLogFlops(2*N-1);
253a40ed3dSBarry Smith   PetscFunctionReturn(0);
261987afe7SBarry Smith }
271987afe7SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense"
308f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
31289bc588SBarry Smith {
324e220ebcSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
33273d9f13SBarry Smith   int          i,N = A->m*A->n,count = 0;
3487828ca2SBarry Smith   PetscScalar  *v = a->v;
353a40ed3dSBarry Smith 
363a40ed3dSBarry Smith   PetscFunctionBegin;
37289bc588SBarry Smith   for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;}
384e220ebcSLois Curfman McInnes 
39273d9f13SBarry Smith   info->rows_global       = (double)A->m;
40273d9f13SBarry Smith   info->columns_global    = (double)A->n;
41273d9f13SBarry Smith   info->rows_local        = (double)A->m;
42273d9f13SBarry Smith   info->columns_local     = (double)A->n;
434e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
444e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
454e220ebcSLois Curfman McInnes   info->nz_used           = (double)count;
464e220ebcSLois Curfman McInnes   info->nz_unneeded       = (double)(N-count);
474e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
484e220ebcSLois Curfman McInnes   info->mallocs           = 0;
494e220ebcSLois Curfman McInnes   info->memory            = A->mem;
504e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
514e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
524e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
534e220ebcSLois Curfman McInnes 
543a40ed3dSBarry Smith   PetscFunctionReturn(0);
55289bc588SBarry Smith }
56289bc588SBarry Smith 
574a2ae208SSatish Balay #undef __FUNCT__
584a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense"
59268466fbSBarry Smith int MatScale_SeqDense(const PetscScalar *alpha,Mat A)
6080cd9d93SLois Curfman McInnes {
61273d9f13SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
62a5ce6ee0Svictorle   int          one = 1,lda = a->lda,j,nz;
6380cd9d93SLois Curfman McInnes 
643a40ed3dSBarry Smith   PetscFunctionBegin;
65a5ce6ee0Svictorle   if (lda>A->m) {
66a5ce6ee0Svictorle     nz = A->m;
67a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
68268466fbSBarry Smith       BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one);
69a5ce6ee0Svictorle     }
70a5ce6ee0Svictorle   } else {
71273d9f13SBarry Smith     nz = A->m*A->n;
72268466fbSBarry Smith     BLscal_(&nz,(PetscScalar*)alpha,a->v,&one);
73a5ce6ee0Svictorle   }
74b0a32e0cSBarry Smith   PetscLogFlops(nz);
753a40ed3dSBarry Smith   PetscFunctionReturn(0);
7680cd9d93SLois Curfman McInnes }
7780cd9d93SLois Curfman McInnes 
78289bc588SBarry Smith /* ---------------------------------------------------------------*/
796831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
806831982aSBarry Smith    rather than put it in the Mat->row slot.*/
814a2ae208SSatish Balay #undef __FUNCT__
824a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense"
83b380c88cSHong Zhang int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo)
84289bc588SBarry Smith {
85a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
8681f479a6Svictorle   PetscFunctionBegin;
87a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
88a07ab388SBarry Smith #else
89c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
90b0a32e0cSBarry Smith   int          info,ierr;
9155659b69SBarry Smith 
923a40ed3dSBarry Smith   PetscFunctionBegin;
93289bc588SBarry Smith   if (!mat->pivots) {
94b0a32e0cSBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr);
95b0a32e0cSBarry Smith     PetscLogObjectMemory(A,A->m*sizeof(int));
96289bc588SBarry Smith   }
97f1af5d2fSBarry Smith   A->factor = FACTOR_LU;
98273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
991b807ce4Svictorle   LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info);
10029bbc08cSBarry Smith   if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization");
10129bbc08cSBarry Smith   if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
102b0a32e0cSBarry Smith   PetscLogFlops((2*A->n*A->n*A->n)/3);
103a07ab388SBarry Smith #endif
1043a40ed3dSBarry Smith   PetscFunctionReturn(0);
105289bc588SBarry Smith }
1066ee01492SSatish Balay 
1074a2ae208SSatish Balay #undef __FUNCT__
1084a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense"
1095609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
11002cad45dSBarry Smith {
11102cad45dSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l;
112a5ce6ee0Svictorle   int          lda = mat->lda,j,m,ierr;
11302cad45dSBarry Smith   Mat          newi;
11402cad45dSBarry Smith 
1153a40ed3dSBarry Smith   PetscFunctionBegin;
1165c5985e7SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&newi);CHKERRQ(ierr);
1175c5985e7SKris Buschelman   ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr);
1185c5985e7SKris Buschelman   ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr);
1195609ef8eSBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
120a5ce6ee0Svictorle     l = (Mat_SeqDense*)newi->data;
121a5ce6ee0Svictorle     if (lda>A->m) {
122a5ce6ee0Svictorle       m = A->m;
123a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
124a5ce6ee0Svictorle 	ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
125a5ce6ee0Svictorle       }
126a5ce6ee0Svictorle     } else {
12787828ca2SBarry Smith       ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
12802cad45dSBarry Smith     }
129a5ce6ee0Svictorle   }
13002cad45dSBarry Smith   newi->assembled = PETSC_TRUE;
13102cad45dSBarry Smith   *newmat = newi;
1323a40ed3dSBarry Smith   PetscFunctionReturn(0);
13302cad45dSBarry Smith }
13402cad45dSBarry Smith 
1354a2ae208SSatish Balay #undef __FUNCT__
1364a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense"
137b380c88cSHong Zhang int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact)
138289bc588SBarry Smith {
1393a40ed3dSBarry Smith   int ierr;
1403a40ed3dSBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
1425609ef8eSBarry Smith   ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
144289bc588SBarry Smith }
1456ee01492SSatish Balay 
1464a2ae208SSatish Balay #undef __FUNCT__
1474a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense"
1488f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact)
149289bc588SBarry Smith {
15002cad45dSBarry Smith   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data;
1511b807ce4Svictorle   int           lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr;
1524482741eSBarry Smith   MatFactorInfo info;
1533a40ed3dSBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
15502cad45dSBarry Smith   /* copy the numerical values */
1561b807ce4Svictorle   if (lda1>m || lda2>m ) {
1571b807ce4Svictorle     for (j=0; j<n; j++) {
1581b807ce4Svictorle       ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1591b807ce4Svictorle     }
1601b807ce4Svictorle   } else {
16187828ca2SBarry Smith     ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1621b807ce4Svictorle   }
16302cad45dSBarry Smith   (*fact)->factor = 0;
1644482741eSBarry Smith   ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr);
1653a40ed3dSBarry Smith   PetscFunctionReturn(0);
166289bc588SBarry Smith }
1676ee01492SSatish Balay 
1684a2ae208SSatish Balay #undef __FUNCT__
1694a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense"
17015e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact)
171289bc588SBarry Smith {
1723a40ed3dSBarry Smith   int ierr;
1733a40ed3dSBarry Smith 
1743a40ed3dSBarry Smith   PetscFunctionBegin;
1753a40ed3dSBarry Smith   ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr);
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177289bc588SBarry Smith }
1786ee01492SSatish Balay 
1794a2ae208SSatish Balay #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense"
18115e8a5b3SHong Zhang int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo)
182289bc588SBarry Smith {
183a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
184a07ab388SBarry Smith   PetscFunctionBegin;
185a07ab388SBarry Smith   SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
186a07ab388SBarry Smith #else
187c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
188606d414cSSatish Balay   int           info,ierr;
18955659b69SBarry Smith 
1903a40ed3dSBarry Smith   PetscFunctionBegin;
191752f5567SLois Curfman McInnes   if (mat->pivots) {
192606d414cSSatish Balay     ierr = PetscFree(mat->pivots);CHKERRQ(ierr);
193b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-A->m*sizeof(int));
194752f5567SLois Curfman McInnes     mat->pivots = 0;
195752f5567SLois Curfman McInnes   }
196273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
1971b807ce4Svictorle   LApotrf_("L",&A->n,mat->v,&mat->lda,&info);
19829bbc08cSBarry Smith   if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1);
199c0bbcb79SLois Curfman McInnes   A->factor = FACTOR_CHOLESKY;
200b0a32e0cSBarry Smith   PetscLogFlops((A->n*A->n*A->n)/3);
201a07ab388SBarry Smith #endif
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203289bc588SBarry Smith }
204289bc588SBarry Smith 
2054a2ae208SSatish Balay #undef __FUNCT__
2060b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense"
2070b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact)
2080b4b3355SBarry Smith {
2090b4b3355SBarry Smith   int ierr;
21015e8a5b3SHong Zhang   MatFactorInfo info;
2110b4b3355SBarry Smith 
2120b4b3355SBarry Smith   PetscFunctionBegin;
21315e8a5b3SHong Zhang   info.fill = 1.0;
21415e8a5b3SHong Zhang   ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr);
2150b4b3355SBarry Smith   PetscFunctionReturn(0);
2160b4b3355SBarry Smith }
2170b4b3355SBarry Smith 
2180b4b3355SBarry Smith #undef __FUNCT__
2194a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense"
2208f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
221289bc588SBarry Smith {
222c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
223a57ad546SLois Curfman McInnes   int          one = 1,info,ierr;
22487828ca2SBarry Smith   PetscScalar  *x,*y;
22567e560aaSBarry Smith 
2263a40ed3dSBarry Smith   PetscFunctionBegin;
227273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2291ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
23087828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
231c0bbcb79SLois Curfman McInnes   if (A->factor == FACTOR_LU) {
232ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
23380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
234ae7cfcebSSatish Balay #else
2351b807ce4Svictorle     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
2362ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
237ae7cfcebSSatish Balay #endif
2387a97a34bSBarry Smith   } else if (A->factor == FACTOR_CHOLESKY){
239ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
24080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
241ae7cfcebSSatish Balay #else
2421b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
2432ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve");
244ae7cfcebSSatish Balay #endif
245289bc588SBarry Smith   }
24629bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
2471ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2481ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
249b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
251289bc588SBarry Smith }
2526ee01492SSatish Balay 
2534a2ae208SSatish Balay #undef __FUNCT__
2544a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense"
2557c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
256da3a660dSBarry Smith {
257c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2587a97a34bSBarry Smith   int          ierr,one = 1,info;
25987828ca2SBarry Smith   PetscScalar  *x,*y;
26067e560aaSBarry Smith 
2613a40ed3dSBarry Smith   PetscFunctionBegin;
262273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
2631ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2641ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
26587828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
266752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
267da3a660dSBarry Smith   if (mat->pivots) {
268ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
26980444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
270ae7cfcebSSatish Balay #else
2711b807ce4Svictorle     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
2722ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
273ae7cfcebSSatish Balay #endif
2747a97a34bSBarry Smith   } else {
275ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
27680444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
277ae7cfcebSSatish Balay #else
2781b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
2792ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
280ae7cfcebSSatish Balay #endif
281da3a660dSBarry Smith   }
2821ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2831ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
284b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n - A->n);
2853a40ed3dSBarry Smith   PetscFunctionReturn(0);
286da3a660dSBarry Smith }
2876ee01492SSatish Balay 
2884a2ae208SSatish Balay #undef __FUNCT__
2894a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense"
2908f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
291da3a660dSBarry Smith {
292c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
2937a97a34bSBarry Smith   int          ierr,one = 1,info;
29487828ca2SBarry Smith   PetscScalar  *x,*y,sone = 1.0;
295da3a660dSBarry Smith   Vec          tmp = 0;
29667e560aaSBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2991ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
300273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
301da3a660dSBarry Smith   if (yy == zz) {
30278b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
303b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
30478b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
305da3a660dSBarry Smith   }
30687828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
307752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
308da3a660dSBarry Smith   if (mat->pivots) {
309ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
31080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
311ae7cfcebSSatish Balay #else
3121b807ce4Svictorle     LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
3132ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
314ae7cfcebSSatish Balay #endif
315a8c6a408SBarry Smith   } else {
316ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
31780444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
318ae7cfcebSSatish Balay #else
3191b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
3202ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
321ae7cfcebSSatish Balay #endif
322da3a660dSBarry Smith   }
323a8c6a408SBarry Smith   if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);}
324a8c6a408SBarry Smith   else     {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);}
3251ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3261ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
327b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
329da3a660dSBarry Smith }
33067e560aaSBarry Smith 
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense"
3337c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
334da3a660dSBarry Smith {
335c0bbcb79SLois Curfman McInnes   Mat_SeqDense  *mat = (Mat_SeqDense*)A->data;
3366abc6512SBarry Smith   int           one = 1,info,ierr;
33787828ca2SBarry Smith   PetscScalar   *x,*y,sone = 1.0;
338da3a660dSBarry Smith   Vec           tmp;
33967e560aaSBarry Smith 
3403a40ed3dSBarry Smith   PetscFunctionBegin;
341273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
3421ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3431ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
344da3a660dSBarry Smith   if (yy == zz) {
34578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
346b0a32e0cSBarry Smith     PetscLogObjectParent(A,tmp);
34778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
348da3a660dSBarry Smith   }
34987828ca2SBarry Smith   ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr);
350752f5567SLois Curfman McInnes   /* assume if pivots exist then use LU; else Cholesky */
351da3a660dSBarry Smith   if (mat->pivots) {
352ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
35380444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
354ae7cfcebSSatish Balay #else
3551b807ce4Svictorle     LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info);
3562ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
357ae7cfcebSSatish Balay #endif
3583a40ed3dSBarry Smith   } else {
359ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
36080444007SBarry Smith     SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
361ae7cfcebSSatish Balay #else
3621b807ce4Svictorle     LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info);
3632ffef20aSBarry Smith     if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve");
364ae7cfcebSSatish Balay #endif
365da3a660dSBarry Smith   }
36690f02eecSBarry Smith   if (tmp) {
36790f02eecSBarry Smith     ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr);
36890f02eecSBarry Smith     ierr = VecDestroy(tmp);CHKERRQ(ierr);
3693a40ed3dSBarry Smith   } else {
37090f02eecSBarry Smith     ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);
37190f02eecSBarry Smith   }
3721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3731ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
374b0a32e0cSBarry Smith   PetscLogFlops(2*A->n*A->n);
3753a40ed3dSBarry Smith   PetscFunctionReturn(0);
376da3a660dSBarry Smith }
377289bc588SBarry Smith /* ------------------------------------------------------------------*/
3784a2ae208SSatish Balay #undef __FUNCT__
3794a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense"
380329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,
381c14dc6b6SHong Zhang                           PetscReal shift,int its,int lits,Vec xx)
382289bc588SBarry Smith {
383c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
38487828ca2SBarry Smith   PetscScalar  *x,*b,*v = mat->v,zero = 0.0,xt;
385273d9f13SBarry Smith   int          ierr,m = A->m,i;
386aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
387bc1b551cSSatish Balay   int          o = 1;
388bc1b551cSSatish Balay #endif
389289bc588SBarry Smith 
3903a40ed3dSBarry Smith   PetscFunctionBegin;
391289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
3923a40ed3dSBarry Smith     /* this is a hack fix, should have another version without the second BLdot */
393bbb6d6a8SBarry Smith     ierr = VecSet(&zero,xx);CHKERRQ(ierr);
394289bc588SBarry Smith   }
3951ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3961ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
397b965ef7fSBarry Smith   its  = its*lits;
39891723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
399289bc588SBarry Smith   while (its--) {
400289bc588SBarry Smith     if (flag & SOR_FORWARD_SWEEP){
401289bc588SBarry Smith       for (i=0; i<m; i++) {
402aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
403f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
404f1747703SBarry Smith            not happy about returning a double complex */
405f1747703SBarry Smith         int         _i;
40687828ca2SBarry Smith         PetscScalar sum = b[i];
407f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4083f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
409f1747703SBarry Smith         }
410f1747703SBarry Smith         xt = sum;
411f1747703SBarry Smith #else
412289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
413f1747703SBarry Smith #endif
41455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
415289bc588SBarry Smith       }
416289bc588SBarry Smith     }
417289bc588SBarry Smith     if (flag & SOR_BACKWARD_SWEEP) {
418289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
419aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
420f1747703SBarry Smith         /* cannot use BLAS dot for complex because compiler/linker is
421f1747703SBarry Smith            not happy about returning a double complex */
422f1747703SBarry Smith         int         _i;
42387828ca2SBarry Smith         PetscScalar sum = b[i];
424f1747703SBarry Smith         for (_i=0; _i<m; _i++) {
4253f6de6efSSatish Balay           sum -= PetscConj(v[i+_i*m])*x[_i];
426f1747703SBarry Smith         }
427f1747703SBarry Smith         xt = sum;
428f1747703SBarry Smith #else
429289bc588SBarry Smith         xt = b[i]-BLdot_(&m,v+i,&m,x,&o);
430f1747703SBarry Smith #endif
43155a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
432289bc588SBarry Smith       }
433289bc588SBarry Smith     }
434289bc588SBarry Smith   }
4351ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
4361ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4373a40ed3dSBarry Smith   PetscFunctionReturn(0);
438289bc588SBarry Smith }
439289bc588SBarry Smith 
440289bc588SBarry Smith /* -----------------------------------------------------------------*/
4414a2ae208SSatish Balay #undef __FUNCT__
4424a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense"
4437c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
444289bc588SBarry Smith {
445c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
44687828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
447ea709b57SSatish Balay   int          ierr,_One=1;
448ea709b57SSatish Balay   PetscScalar  _DOne=1.0,_DZero=0.0;
4493a40ed3dSBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
451273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4531ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
4541b807ce4Svictorle   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4551ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4561ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
457b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->n);
4583a40ed3dSBarry Smith   PetscFunctionReturn(0);
459289bc588SBarry Smith }
4606ee01492SSatish Balay 
4614a2ae208SSatish Balay #undef __FUNCT__
4624a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense"
46344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy)
464289bc588SBarry Smith {
465c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
46687828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0;
467329f5518SBarry Smith   int          ierr,_One=1;
4683a40ed3dSBarry Smith 
4693a40ed3dSBarry Smith   PetscFunctionBegin;
470273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
4711ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4721ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
4731b807ce4Svictorle   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One);
4741ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4751ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
476b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n - A->m);
4773a40ed3dSBarry Smith   PetscFunctionReturn(0);
478289bc588SBarry Smith }
4796ee01492SSatish Balay 
4804a2ae208SSatish Balay #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense"
48244cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
483289bc588SBarry Smith {
484c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
48587828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y,_DOne=1.0;
486329f5518SBarry Smith   int          ierr,_One=1;
4873a40ed3dSBarry Smith 
4883a40ed3dSBarry Smith   PetscFunctionBegin;
489273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
490600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
4911ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4921ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
4931b807ce4Svictorle   LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
4941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4951ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
496b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
4973a40ed3dSBarry Smith   PetscFunctionReturn(0);
498289bc588SBarry Smith }
4996ee01492SSatish Balay 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense"
5027c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
503289bc588SBarry Smith {
504c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
50587828ca2SBarry Smith   PetscScalar  *v = mat->v,*x,*y;
5067a97a34bSBarry Smith   int          ierr,_One=1;
50787828ca2SBarry Smith   PetscScalar  _DOne=1.0;
5083a40ed3dSBarry Smith 
5093a40ed3dSBarry Smith   PetscFunctionBegin;
510273d9f13SBarry Smith   if (!A->m || !A->n) PetscFunctionReturn(0);
511600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
5121ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5131ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5141b807ce4Svictorle   LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One);
5151ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5161ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
517b0a32e0cSBarry Smith   PetscLogFlops(2*A->m*A->n);
5183a40ed3dSBarry Smith   PetscFunctionReturn(0);
519289bc588SBarry Smith }
520289bc588SBarry Smith 
521289bc588SBarry Smith /* -----------------------------------------------------------------*/
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense"
52487828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
525289bc588SBarry Smith {
526c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
52787828ca2SBarry Smith   PetscScalar  *v;
528b0a32e0cSBarry Smith   int          i,ierr;
52967e560aaSBarry Smith 
5303a40ed3dSBarry Smith   PetscFunctionBegin;
531273d9f13SBarry Smith   *ncols = A->n;
532289bc588SBarry Smith   if (cols) {
533b0a32e0cSBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr);
534273d9f13SBarry Smith     for (i=0; i<A->n; i++) (*cols)[i] = i;
535289bc588SBarry Smith   }
536289bc588SBarry Smith   if (vals) {
53787828ca2SBarry Smith     ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr);
538289bc588SBarry Smith     v    = mat->v + row;
5391b807ce4Svictorle     for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;}
540289bc588SBarry Smith   }
5413a40ed3dSBarry Smith   PetscFunctionReturn(0);
542289bc588SBarry Smith }
5436ee01492SSatish Balay 
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense"
54687828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
547289bc588SBarry Smith {
548606d414cSSatish Balay   int ierr;
549606d414cSSatish Balay   PetscFunctionBegin;
550606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
551606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
5523a40ed3dSBarry Smith   PetscFunctionReturn(0);
553289bc588SBarry Smith }
554289bc588SBarry Smith /* ----------------------------------------------------------------*/
5554a2ae208SSatish Balay #undef __FUNCT__
5564a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense"
557f15d580aSBarry Smith int MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv)
558289bc588SBarry Smith {
559c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
560cddbea37SSatish Balay   int          i,j,idx=0;
561d6dfbf8fSBarry Smith 
5623a40ed3dSBarry Smith   PetscFunctionBegin;
563289bc588SBarry Smith   if (!mat->roworiented) {
564dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
565289bc588SBarry Smith       for (j=0; j<n; j++) {
566cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
567aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
568590ac198SBarry Smith         if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1);
56958804f6dSBarry Smith #endif
570289bc588SBarry Smith         for (i=0; i<m; i++) {
571cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
572aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5735c8e7597SSatish Balay           if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1);
57458804f6dSBarry Smith #endif
575cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
576289bc588SBarry Smith         }
577289bc588SBarry Smith       }
5783a40ed3dSBarry Smith     } else {
579289bc588SBarry Smith       for (j=0; j<n; j++) {
580cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
581aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
582590ac198SBarry 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;}
586aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5875c8e7597SSatish Balay           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       }
592289bc588SBarry Smith     }
5933a40ed3dSBarry Smith   } else {
594dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
595e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
596cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
597aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
5985c8e7597SSatish Balay         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
59958804f6dSBarry Smith #endif
600e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
601cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
602aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6035c8e7597SSatish Balay           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
60458804f6dSBarry Smith #endif
605cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
606e8d4e0b9SBarry Smith         }
607e8d4e0b9SBarry Smith       }
6083a40ed3dSBarry Smith     } else {
609289bc588SBarry Smith       for (i=0; i<m; i++) {
610cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
611aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6125c8e7597SSatish Balay         if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1);
61358804f6dSBarry Smith #endif
614289bc588SBarry Smith         for (j=0; j<n; j++) {
615cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
616aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
6175c8e7597SSatish Balay           if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1);
61858804f6dSBarry Smith #endif
619cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
620289bc588SBarry Smith         }
621289bc588SBarry Smith       }
622289bc588SBarry Smith     }
623e8d4e0b9SBarry Smith   }
6243a40ed3dSBarry Smith   PetscFunctionReturn(0);
625289bc588SBarry Smith }
626e8d4e0b9SBarry Smith 
6274a2ae208SSatish Balay #undef __FUNCT__
6284a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense"
629f15d580aSBarry Smith int MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[])
630ae80bb75SLois Curfman McInnes {
631ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
632ae80bb75SLois Curfman McInnes   int          i,j;
63387828ca2SBarry Smith   PetscScalar  *vpt = v;
634ae80bb75SLois Curfman McInnes 
6353a40ed3dSBarry Smith   PetscFunctionBegin;
636ae80bb75SLois Curfman McInnes   /* row-oriented output */
637ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
638ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
6391b807ce4Svictorle       *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]];
640ae80bb75SLois Curfman McInnes     }
641ae80bb75SLois Curfman McInnes   }
6423a40ed3dSBarry Smith   PetscFunctionReturn(0);
643ae80bb75SLois Curfman McInnes }
644ae80bb75SLois Curfman McInnes 
645289bc588SBarry Smith /* -----------------------------------------------------------------*/
646289bc588SBarry Smith 
647e090d566SSatish Balay #include "petscsys.h"
64888e32edaSLois Curfman McInnes 
6494a2ae208SSatish Balay #undef __FUNCT__
6504a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense"
6518e9aea5cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A)
65288e32edaSLois Curfman McInnes {
65388e32edaSLois Curfman McInnes   Mat_SeqDense *a;
65488e32edaSLois Curfman McInnes   Mat          B;
65588e32edaSLois Curfman McInnes   int          *scols,i,j,nz,ierr,fd,header[4],size;
65688e32edaSLois Curfman McInnes   int          *rowlengths = 0,M,N,*cols;
65787828ca2SBarry Smith   PetscScalar  *vals,*svals,*v,*w;
65819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
65988e32edaSLois Curfman McInnes 
6603a40ed3dSBarry Smith   PetscFunctionBegin;
661d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
66229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
663b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
6640752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
665552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
66688e32edaSLois Curfman McInnes   M = header[1]; N = header[2]; nz = header[3];
66788e32edaSLois Curfman McInnes 
66890ace30eSBarry Smith   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
6695c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
670be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6715c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
67290ace30eSBarry Smith     B    = *A;
67390ace30eSBarry Smith     a    = (Mat_SeqDense*)B->data;
6744a41a4d3SSatish Balay     v    = a->v;
6754a41a4d3SSatish Balay     /* Allocate some temp space to read in the values and then flip them
6764a41a4d3SSatish Balay        from row major to column major */
677b72907f3SBarry Smith     ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr);
67890ace30eSBarry Smith     /* read in nonzero values */
6794a41a4d3SSatish Balay     ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr);
6804a41a4d3SSatish Balay     /* now flip the values and store them in the matrix*/
6814a41a4d3SSatish Balay     for (j=0; j<N; j++) {
6824a41a4d3SSatish Balay       for (i=0; i<M; i++) {
6834a41a4d3SSatish Balay         *v++ =w[i*N+j];
6844a41a4d3SSatish Balay       }
6854a41a4d3SSatish Balay     }
686606d414cSSatish Balay     ierr = PetscFree(w);CHKERRQ(ierr);
6876d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6886d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68990ace30eSBarry Smith   } else {
69088e32edaSLois Curfman McInnes     /* read row lengths */
691b0a32e0cSBarry Smith     ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr);
6920752156aSBarry Smith     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
69388e32edaSLois Curfman McInnes 
69488e32edaSLois Curfman McInnes     /* create our matrix */
6955c5985e7SKris Buschelman     ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr);
696be5d1d56SKris Buschelman     ierr = MatSetType(*A,type);CHKERRQ(ierr);
6975c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr);
69888e32edaSLois Curfman McInnes     B = *A;
69988e32edaSLois Curfman McInnes     a = (Mat_SeqDense*)B->data;
70088e32edaSLois Curfman McInnes     v = a->v;
70188e32edaSLois Curfman McInnes 
70288e32edaSLois Curfman McInnes     /* read column indices and nonzeros */
703b0a32e0cSBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr);
704b0a32e0cSBarry Smith     cols = scols;
7050752156aSBarry Smith     ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
70687828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr);
707b0a32e0cSBarry Smith     vals = svals;
7080752156aSBarry Smith     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
70988e32edaSLois Curfman McInnes 
71088e32edaSLois Curfman McInnes     /* insert into matrix */
71188e32edaSLois Curfman McInnes     for (i=0; i<M; i++) {
71288e32edaSLois Curfman McInnes       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
71388e32edaSLois Curfman McInnes       svals += rowlengths[i]; scols += rowlengths[i];
71488e32edaSLois Curfman McInnes     }
715606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
716606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
717606d414cSSatish Balay     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
71888e32edaSLois Curfman McInnes 
7196d4a8577SBarry Smith     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7206d4a8577SBarry Smith     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72190ace30eSBarry Smith   }
7223a40ed3dSBarry Smith   PetscFunctionReturn(0);
72388e32edaSLois Curfman McInnes }
72488e32edaSLois Curfman McInnes 
725e090d566SSatish Balay #include "petscsys.h"
726932b0c3eSLois Curfman McInnes 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII"
729b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
730289bc588SBarry Smith {
731932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
732fb9695e5SSatish Balay   int               ierr,i,j;
733fb9695e5SSatish Balay   char              *name;
73487828ca2SBarry Smith   PetscScalar       *v;
735f3ef73ceSBarry Smith   PetscViewerFormat format;
736932b0c3eSLois Curfman McInnes 
7373a40ed3dSBarry Smith   PetscFunctionBegin;
738435da068SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
739b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
740456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
7413a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
742fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
743b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
744273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
74544cd7ae7SLois Curfman McInnes       v = a->v + i;
746b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr);
747273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
748aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
749329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
75062b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
751329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
75262b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr);
7536831982aSBarry Smith         }
75480cd9d93SLois Curfman McInnes #else
7556831982aSBarry Smith         if (*v) {
75662b941d6SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr);
7576831982aSBarry Smith         }
75880cd9d93SLois Curfman McInnes #endif
7591b807ce4Svictorle         v += a->lda;
76080cd9d93SLois Curfman McInnes       }
761b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
76280cd9d93SLois Curfman McInnes     }
763b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
7643a40ed3dSBarry Smith   } else {
765b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
766aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
767ffac6cdbSBarry Smith     PetscTruth allreal = PETSC_TRUE;
76847989497SBarry Smith     /* determine if matrix has all real values */
76947989497SBarry Smith     v = a->v;
770201f5f94SSatish Balay     for (i=0; i<A->m*A->n; i++) {
771ffac6cdbSBarry Smith 	if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;}
77247989497SBarry Smith     }
77347989497SBarry Smith #endif
774fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
7753a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
776b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr);
777fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr);
778fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
779ffac6cdbSBarry Smith     }
780ffac6cdbSBarry Smith 
781273d9f13SBarry Smith     for (i=0; i<A->m; i++) {
782932b0c3eSLois Curfman McInnes       v = a->v + i;
783273d9f13SBarry Smith       for (j=0; j<A->n; j++) {
784aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
78547989497SBarry Smith         if (allreal) {
7861b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr);
78747989497SBarry Smith         } else {
7881b807ce4Svictorle           ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr);
78947989497SBarry Smith         }
790289bc588SBarry Smith #else
7911b807ce4Svictorle         ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr);
792289bc588SBarry Smith #endif
7931b807ce4Svictorle 	v += a->lda;
794289bc588SBarry Smith       }
795b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
796289bc588SBarry Smith     }
797fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
798b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
799ffac6cdbSBarry Smith     }
800b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
801da3a660dSBarry Smith   }
802b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
804289bc588SBarry Smith }
805289bc588SBarry Smith 
8064a2ae208SSatish Balay #undef __FUNCT__
8074a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary"
808b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
809932b0c3eSLois Curfman McInnes {
810932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
811273d9f13SBarry Smith   int               ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n;
81287828ca2SBarry Smith   PetscScalar       *v,*anonz,*vals;
813f3ef73ceSBarry Smith   PetscViewerFormat format;
814932b0c3eSLois Curfman McInnes 
8153a40ed3dSBarry Smith   PetscFunctionBegin;
816b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
81790ace30eSBarry Smith 
818b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
819fb9695e5SSatish Balay   if (format == PETSC_VIEWER_BINARY_NATIVE) {
82090ace30eSBarry Smith     /* store the matrix as a dense matrix */
821b0a32e0cSBarry Smith     ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr);
822552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
82390ace30eSBarry Smith     col_lens[1] = m;
82490ace30eSBarry Smith     col_lens[2] = n;
82590ace30eSBarry Smith     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
8260752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr);
827606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
82890ace30eSBarry Smith 
82990ace30eSBarry Smith     /* write out matrix, by rows */
83087828ca2SBarry Smith     ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
83190ace30eSBarry Smith     v    = a->v;
83290ace30eSBarry Smith     for (i=0; i<m; i++) {
83390ace30eSBarry Smith       for (j=0; j<n; j++) {
83490ace30eSBarry Smith         vals[i + j*m] = *v++;
83590ace30eSBarry Smith       }
83690ace30eSBarry Smith     }
8370752156aSBarry Smith     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr);
838606d414cSSatish Balay     ierr = PetscFree(vals);CHKERRQ(ierr);
83990ace30eSBarry Smith   } else {
840b0a32e0cSBarry Smith     ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr);
841552e946dSBarry Smith     col_lens[0] = MAT_FILE_COOKIE;
842932b0c3eSLois Curfman McInnes     col_lens[1] = m;
843932b0c3eSLois Curfman McInnes     col_lens[2] = n;
844932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
845932b0c3eSLois Curfman McInnes 
846932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
847932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
8480752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr);
849932b0c3eSLois Curfman McInnes 
850932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
851932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
852932b0c3eSLois Curfman McInnes     ict = 0;
853932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
854932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
855932b0c3eSLois Curfman McInnes     }
8560752156aSBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr);
857606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
858932b0c3eSLois Curfman McInnes 
859932b0c3eSLois Curfman McInnes     /* store nonzero values */
86087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr);
861932b0c3eSLois Curfman McInnes     ict  = 0;
862932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
863932b0c3eSLois Curfman McInnes       v = a->v + i;
864932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
8651b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
866932b0c3eSLois Curfman McInnes       }
867932b0c3eSLois Curfman McInnes     }
8680752156aSBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr);
869606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
87090ace30eSBarry Smith   }
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
872932b0c3eSLois Curfman McInnes }
873932b0c3eSLois Curfman McInnes 
8744a2ae208SSatish Balay #undef __FUNCT__
8754a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom"
876b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
877f1af5d2fSBarry Smith {
878f1af5d2fSBarry Smith   Mat               A = (Mat) Aa;
879f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
880fb9695e5SSatish Balay   int               m = A->m,n = A->n,color,i,j,ierr;
88187828ca2SBarry Smith   PetscScalar       *v = a->v;
882b0a32e0cSBarry Smith   PetscViewer       viewer;
883b0a32e0cSBarry Smith   PetscDraw         popup;
884329f5518SBarry Smith   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0;
885f3ef73ceSBarry Smith   PetscViewerFormat format;
886f1af5d2fSBarry Smith 
887f1af5d2fSBarry Smith   PetscFunctionBegin;
888f1af5d2fSBarry Smith 
889f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
890b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
891b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
892f1af5d2fSBarry Smith 
893f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
894fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
895f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
896b0a32e0cSBarry Smith     color = PETSC_DRAW_BLUE;
897f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
898f1af5d2fSBarry Smith       x_l = j;
899f1af5d2fSBarry Smith       x_r = x_l + 1.0;
900f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
901f1af5d2fSBarry Smith         y_l = m - i - 1.0;
902f1af5d2fSBarry Smith         y_r = y_l + 1.0;
903f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX)
904329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
905b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
906329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
907b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
908f1af5d2fSBarry Smith         } else {
909f1af5d2fSBarry Smith           continue;
910f1af5d2fSBarry Smith         }
911f1af5d2fSBarry Smith #else
912f1af5d2fSBarry Smith         if (v[j*m+i] >  0.) {
913b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
914f1af5d2fSBarry Smith         } else if (v[j*m+i] <  0.) {
915b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
916f1af5d2fSBarry Smith         } else {
917f1af5d2fSBarry Smith           continue;
918f1af5d2fSBarry Smith         }
919f1af5d2fSBarry Smith #endif
920b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
921f1af5d2fSBarry Smith       }
922f1af5d2fSBarry Smith     }
923f1af5d2fSBarry Smith   } else {
924f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
925f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
926f1af5d2fSBarry Smith     for(i = 0; i < m*n; i++) {
927f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
928f1af5d2fSBarry Smith     }
929b0a32e0cSBarry Smith     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
930b0a32e0cSBarry Smith     ierr  = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
931b0a32e0cSBarry Smith     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);}
932f1af5d2fSBarry Smith     for(j = 0; j < n; j++) {
933f1af5d2fSBarry Smith       x_l = j;
934f1af5d2fSBarry Smith       x_r = x_l + 1.0;
935f1af5d2fSBarry Smith       for(i = 0; i < m; i++) {
936f1af5d2fSBarry Smith         y_l   = m - i - 1.0;
937f1af5d2fSBarry Smith         y_r   = y_l + 1.0;
938b0a32e0cSBarry Smith         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i]));
939b0a32e0cSBarry Smith         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
940f1af5d2fSBarry Smith       }
941f1af5d2fSBarry Smith     }
942f1af5d2fSBarry Smith   }
943f1af5d2fSBarry Smith   PetscFunctionReturn(0);
944f1af5d2fSBarry Smith }
945f1af5d2fSBarry Smith 
9464a2ae208SSatish Balay #undef __FUNCT__
9474a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw"
948b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
949f1af5d2fSBarry Smith {
950b0a32e0cSBarry Smith   PetscDraw  draw;
951f1af5d2fSBarry Smith   PetscTruth isnull;
952329f5518SBarry Smith   PetscReal  xr,yr,xl,yl,h,w;
953f1af5d2fSBarry Smith   int        ierr;
954f1af5d2fSBarry Smith 
955f1af5d2fSBarry Smith   PetscFunctionBegin;
956b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
957b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
958f1af5d2fSBarry Smith   if (isnull == PETSC_TRUE) PetscFunctionReturn(0);
959f1af5d2fSBarry Smith 
960f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
961273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
962f1af5d2fSBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
963b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
964b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
965f1af5d2fSBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
966f1af5d2fSBarry Smith   PetscFunctionReturn(0);
967f1af5d2fSBarry Smith }
968f1af5d2fSBarry Smith 
9694a2ae208SSatish Balay #undef __FUNCT__
9704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense"
971b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer)
972932b0c3eSLois Curfman McInnes {
973932b0c3eSLois Curfman McInnes   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
974bcd2baecSBarry Smith   int          ierr;
97532077d6dSBarry Smith   PetscTruth   issocket,iascii,isbinary,isdraw;
976932b0c3eSLois Curfman McInnes 
9773a40ed3dSBarry Smith   PetscFunctionBegin;
978b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
97932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
980fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
981fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
9820f5bd95cSBarry Smith 
9830f5bd95cSBarry Smith   if (issocket) {
984*634064b4SBarry Smith     if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA");
985b0a32e0cSBarry Smith     ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr);
98632077d6dSBarry Smith   } else if (iascii) {
9873a40ed3dSBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
9880f5bd95cSBarry Smith   } else if (isbinary) {
9893a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
990f1af5d2fSBarry Smith   } else if (isdraw) {
991f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
9925cd90555SBarry Smith   } else {
993958c9bccSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name);
994932b0c3eSLois Curfman McInnes   }
9953a40ed3dSBarry Smith   PetscFunctionReturn(0);
996932b0c3eSLois Curfman McInnes }
997289bc588SBarry Smith 
9984a2ae208SSatish Balay #undef __FUNCT__
9994a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense"
1000c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat)
1001289bc588SBarry Smith {
1002ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)mat->data;
100390f02eecSBarry Smith   int          ierr;
100490f02eecSBarry Smith 
10053a40ed3dSBarry Smith   PetscFunctionBegin;
1006aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1007b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n);
1008a5a9c739SBarry Smith #endif
1009606d414cSSatish Balay   if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);}
1010606d414cSSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1011606d414cSSatish Balay   ierr = PetscFree(l);CHKERRQ(ierr);
10123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1013289bc588SBarry Smith }
1014289bc588SBarry Smith 
10154a2ae208SSatish Balay #undef __FUNCT__
10164a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense"
10178f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout)
1018289bc588SBarry Smith {
1019c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10201b807ce4Svictorle   int          k,j,m,n,M,ierr;
102187828ca2SBarry Smith   PetscScalar  *v,tmp;
102248b35521SBarry Smith 
10233a40ed3dSBarry Smith   PetscFunctionBegin;
10241b807ce4Svictorle   v = mat->v; m = A->m; M = mat->lda; n = A->n;
10257c922b88SBarry Smith   if (!matout) { /* in place transpose */
1026a5ce6ee0Svictorle     if (m != n) {
1027*634064b4SBarry Smith       SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place");
1028d64ed03dSBarry Smith     } else {
1029d3e5ee88SLois Curfman McInnes       for (j=0; j<m; j++) {
1030289bc588SBarry Smith         for (k=0; k<j; k++) {
10311b807ce4Svictorle           tmp = v[j + k*M];
10321b807ce4Svictorle           v[j + k*M] = v[k + j*M];
10331b807ce4Svictorle           v[k + j*M] = tmp;
1034289bc588SBarry Smith         }
1035289bc588SBarry Smith       }
1036d64ed03dSBarry Smith     }
10373a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1038d3e5ee88SLois Curfman McInnes     Mat          tmat;
1039ec8511deSBarry Smith     Mat_SeqDense *tmatd;
104087828ca2SBarry Smith     PetscScalar  *v2;
1041ea709b57SSatish Balay 
10425c5985e7SKris Buschelman     ierr  = MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);CHKERRQ(ierr);
10435c5985e7SKris Buschelman     ierr  = MatSetType(tmat,A->type_name);CHKERRQ(ierr);
10445c5985e7SKris Buschelman     ierr  = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr);
1045ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
10460de55854SLois Curfman McInnes     v = mat->v; v2 = tmatd->v;
1047d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
10481b807ce4Svictorle       for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M];
1049d3e5ee88SLois Curfman McInnes     }
10506d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10516d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1052d3e5ee88SLois Curfman McInnes     *matout = tmat;
105348b35521SBarry Smith   }
10543a40ed3dSBarry Smith   PetscFunctionReturn(0);
1055289bc588SBarry Smith }
1056289bc588SBarry Smith 
10574a2ae208SSatish Balay #undef __FUNCT__
10584a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense"
10598f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg)
1060289bc588SBarry Smith {
1061c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1062c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
1063a43ee2ecSKris Buschelman   int          i,j;
106487828ca2SBarry Smith   PetscScalar  *v1 = mat1->v,*v2 = mat2->v;
10659ea5d5aeSSatish Balay 
10663a40ed3dSBarry Smith   PetscFunctionBegin;
1067273d9f13SBarry Smith   if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1068273d9f13SBarry Smith   if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10691b807ce4Svictorle   for (i=0; i<A1->m; i++) {
10701b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
10711b807ce4Svictorle     for (j=0; j<A1->n; j++) {
10723a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
10731b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
10741b807ce4Svictorle     }
1075289bc588SBarry Smith   }
107677c4ece6SBarry Smith   *flg = PETSC_TRUE;
10773a40ed3dSBarry Smith   PetscFunctionReturn(0);
1078289bc588SBarry Smith }
1079289bc588SBarry Smith 
10804a2ae208SSatish Balay #undef __FUNCT__
10814a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense"
10828f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v)
1083289bc588SBarry Smith {
1084c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
10857a97a34bSBarry Smith   int          ierr,i,n,len;
108687828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
108744cd7ae7SLois Curfman McInnes 
10883a40ed3dSBarry Smith   PetscFunctionBegin;
10897a97a34bSBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
10907a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
10911ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1092273d9f13SBarry Smith   len = PetscMin(A->m,A->n);
1093273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
109444cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
10951b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1096289bc588SBarry Smith   }
10971ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
10983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1099289bc588SBarry Smith }
1100289bc588SBarry Smith 
11014a2ae208SSatish Balay #undef __FUNCT__
11024a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense"
11038f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1104289bc588SBarry Smith {
1105c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
110687828ca2SBarry Smith   PetscScalar  *l,*r,x,*v;
1107273d9f13SBarry Smith   int          ierr,i,j,m = A->m,n = A->n;
110855659b69SBarry Smith 
11093a40ed3dSBarry Smith   PetscFunctionBegin;
111028988994SBarry Smith   if (ll) {
11117a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
11121ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1113273d9f13SBarry Smith     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1114da3a660dSBarry Smith     for (i=0; i<m; i++) {
1115da3a660dSBarry Smith       x = l[i];
1116da3a660dSBarry Smith       v = mat->v + i;
1117da3a660dSBarry Smith       for (j=0; j<n; j++) { (*v) *= x; v+= m;}
1118da3a660dSBarry Smith     }
11191ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1120b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1121da3a660dSBarry Smith   }
112228988994SBarry Smith   if (rr) {
11237a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
11241ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1125273d9f13SBarry Smith     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1126da3a660dSBarry Smith     for (i=0; i<n; i++) {
1127da3a660dSBarry Smith       x = r[i];
1128da3a660dSBarry Smith       v = mat->v + i*m;
1129da3a660dSBarry Smith       for (j=0; j<m; j++) { (*v++) *= x;}
1130da3a660dSBarry Smith     }
11311ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1132b0a32e0cSBarry Smith     PetscLogFlops(n*m);
1133da3a660dSBarry Smith   }
11343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1135289bc588SBarry Smith }
1136289bc588SBarry Smith 
11374a2ae208SSatish Balay #undef __FUNCT__
11384a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense"
1139064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1140289bc588SBarry Smith {
1141c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
114287828ca2SBarry Smith   PetscScalar  *v = mat->v;
1143329f5518SBarry Smith   PetscReal    sum = 0.0;
1144a5ce6ee0Svictorle   int          lda=mat->lda,m=A->m,i,j;
114555659b69SBarry Smith 
11463a40ed3dSBarry Smith   PetscFunctionBegin;
1147289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1148a5ce6ee0Svictorle     if (lda>m) {
1149a5ce6ee0Svictorle       for (j=0; j<A->n; j++) {
1150a5ce6ee0Svictorle 	v = mat->v+j*lda;
1151a5ce6ee0Svictorle 	for (i=0; i<m; i++) {
1152a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX)
1153a5ce6ee0Svictorle 	  sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1154a5ce6ee0Svictorle #else
1155a5ce6ee0Svictorle 	  sum += (*v)*(*v); v++;
1156a5ce6ee0Svictorle #endif
1157a5ce6ee0Svictorle 	}
1158a5ce6ee0Svictorle       }
1159a5ce6ee0Svictorle     } else {
1160273d9f13SBarry Smith       for (i=0; i<A->n*A->m; i++) {
1161aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1162329f5518SBarry Smith 	sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1163289bc588SBarry Smith #else
1164289bc588SBarry Smith 	sum += (*v)*(*v); v++;
1165289bc588SBarry Smith #endif
1166289bc588SBarry Smith       }
1167a5ce6ee0Svictorle     }
1168064f8208SBarry Smith     *nrm = sqrt(sum);
1169b0a32e0cSBarry Smith     PetscLogFlops(2*A->n*A->m);
11703a40ed3dSBarry Smith   } else if (type == NORM_1) {
1171064f8208SBarry Smith     *nrm = 0.0;
1172273d9f13SBarry Smith     for (j=0; j<A->n; j++) {
11731b807ce4Svictorle       v = mat->v + j*mat->lda;
1174289bc588SBarry Smith       sum = 0.0;
1175273d9f13SBarry Smith       for (i=0; i<A->m; i++) {
117633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1177289bc588SBarry Smith       }
1178064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1179289bc588SBarry Smith     }
1180b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11813a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1182064f8208SBarry Smith     *nrm = 0.0;
1183273d9f13SBarry Smith     for (j=0; j<A->m; j++) {
1184289bc588SBarry Smith       v = mat->v + j;
1185289bc588SBarry Smith       sum = 0.0;
1186273d9f13SBarry Smith       for (i=0; i<A->n; i++) {
11871b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1188289bc588SBarry Smith       }
1189064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1190289bc588SBarry Smith     }
1191b0a32e0cSBarry Smith     PetscLogFlops(A->n*A->m);
11923a40ed3dSBarry Smith   } else {
119329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No two norm");
1194289bc588SBarry Smith   }
11953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1196289bc588SBarry Smith }
1197289bc588SBarry Smith 
11984a2ae208SSatish Balay #undef __FUNCT__
11994a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense"
12008f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op)
1201289bc588SBarry Smith {
1202c0bbcb79SLois Curfman McInnes   Mat_SeqDense *aij = (Mat_SeqDense*)A->data;
120367e560aaSBarry Smith 
12043a40ed3dSBarry Smith   PetscFunctionBegin;
1205b5a2b587SKris Buschelman   switch (op) {
1206b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
1207b5a2b587SKris Buschelman     aij->roworiented = PETSC_TRUE;
1208b5a2b587SKris Buschelman     break;
1209b5a2b587SKris Buschelman   case MAT_COLUMN_ORIENTED:
1210b5a2b587SKris Buschelman     aij->roworiented = PETSC_FALSE;
1211b5a2b587SKris Buschelman     break;
1212b5a2b587SKris Buschelman   case MAT_ROWS_SORTED:
1213b5a2b587SKris Buschelman   case MAT_ROWS_UNSORTED:
1214b5a2b587SKris Buschelman   case MAT_COLUMNS_SORTED:
1215b5a2b587SKris Buschelman   case MAT_COLUMNS_UNSORTED:
1216b5a2b587SKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1217b5a2b587SKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1218b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1219b5a2b587SKris Buschelman   case MAT_NO_NEW_DIAGONALS:
1220b5a2b587SKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1221b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1222b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
1223b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n");
1224b5a2b587SKris Buschelman     break;
122577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
122677e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
12279a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
12289a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
12299a4540c5SBarry Smith   case MAT_HERMITIAN:
12309a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
12319a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
12329a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
123377e54ba9SKris Buschelman     break;
1234b5a2b587SKris Buschelman   default:
123529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
12363a40ed3dSBarry Smith   }
12373a40ed3dSBarry Smith   PetscFunctionReturn(0);
1238289bc588SBarry Smith }
1239289bc588SBarry Smith 
12404a2ae208SSatish Balay #undef __FUNCT__
12414a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense"
12428f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A)
12436f0a148fSBarry Smith {
1244ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1245a5ce6ee0Svictorle   int          lda=l->lda,m=A->m,j, ierr;
12463a40ed3dSBarry Smith 
12473a40ed3dSBarry Smith   PetscFunctionBegin;
1248a5ce6ee0Svictorle   if (lda>m) {
1249a5ce6ee0Svictorle     for (j=0; j<A->n; j++) {
1250a5ce6ee0Svictorle       ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr);
1251a5ce6ee0Svictorle     }
1252a5ce6ee0Svictorle   } else {
125387828ca2SBarry Smith     ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1254a5ce6ee0Svictorle   }
12553a40ed3dSBarry Smith   PetscFunctionReturn(0);
12566f0a148fSBarry Smith }
12576f0a148fSBarry Smith 
12584a2ae208SSatish Balay #undef __FUNCT__
12594a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense"
12608f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs)
12614e220ebcSLois Curfman McInnes {
12623a40ed3dSBarry Smith   PetscFunctionBegin;
12634e220ebcSLois Curfman McInnes   *bs = 1;
12643a40ed3dSBarry Smith   PetscFunctionReturn(0);
12654e220ebcSLois Curfman McInnes }
12664e220ebcSLois Curfman McInnes 
12674a2ae208SSatish Balay #undef __FUNCT__
12684a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense"
1269268466fbSBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag)
12706f0a148fSBarry Smith {
1271ec8511deSBarry Smith   Mat_SeqDense *l = (Mat_SeqDense*)A->data;
1272273d9f13SBarry Smith   int          n = A->n,i,j,ierr,N,*rows;
127387828ca2SBarry Smith   PetscScalar  *slot;
127455659b69SBarry Smith 
12753a40ed3dSBarry Smith   PetscFunctionBegin;
1276e03a110bSBarry Smith   ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr);
127778b31e54SBarry Smith   ierr = ISGetIndices(is,&rows);CHKERRQ(ierr);
12786f0a148fSBarry Smith   for (i=0; i<N; i++) {
12796f0a148fSBarry Smith     slot = l->v + rows[i];
12806f0a148fSBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += n;}
12816f0a148fSBarry Smith   }
12826f0a148fSBarry Smith   if (diag) {
12836f0a148fSBarry Smith     for (i=0; i<N; i++) {
12846f0a148fSBarry Smith       slot = l->v + (n+1)*rows[i];
1285260925b8SBarry Smith       *slot = *diag;
12866f0a148fSBarry Smith     }
12876f0a148fSBarry Smith   }
1288260925b8SBarry Smith   ISRestoreIndices(is,&rows);
12893a40ed3dSBarry Smith   PetscFunctionReturn(0);
12906f0a148fSBarry Smith }
1291557bce09SLois Curfman McInnes 
12924a2ae208SSatish Balay #undef __FUNCT__
12934a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense"
12944e7234bfSBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar *array[])
129564e87e97SBarry Smith {
1296c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
12973a40ed3dSBarry Smith 
12983a40ed3dSBarry Smith   PetscFunctionBegin;
129964e87e97SBarry Smith   *array = mat->v;
13003a40ed3dSBarry Smith   PetscFunctionReturn(0);
130164e87e97SBarry Smith }
13020754003eSLois Curfman McInnes 
13034a2ae208SSatish Balay #undef __FUNCT__
13044a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense"
13054e7234bfSBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1306ff14e315SSatish Balay {
13073a40ed3dSBarry Smith   PetscFunctionBegin;
130809b544d4SBarry Smith   *array = 0; /* user cannot accidently use the array later */
13093a40ed3dSBarry Smith   PetscFunctionReturn(0);
1310ff14e315SSatish Balay }
13110754003eSLois Curfman McInnes 
13124a2ae208SSatish Balay #undef __FUNCT__
13134a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense"
13147b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
13150754003eSLois Curfman McInnes {
1316c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
1317273d9f13SBarry Smith   int          i,j,ierr,m = A->m,*irow,*icol,nrows,ncols;
131887828ca2SBarry Smith   PetscScalar  *av,*bv,*v = mat->v;
13190754003eSLois Curfman McInnes   Mat          newmat;
13200754003eSLois Curfman McInnes 
13213a40ed3dSBarry Smith   PetscFunctionBegin;
132278b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
132378b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1324e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1325e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
13260754003eSLois Curfman McInnes 
1327182d2002SSatish Balay   /* Check submatrixcall */
1328182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
1329182d2002SSatish Balay     int n_cols,n_rows;
1330182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
133129bbc08cSBarry Smith     if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1332182d2002SSatish Balay     newmat = *B;
1333182d2002SSatish Balay   } else {
13340754003eSLois Curfman McInnes     /* Create and fill new matrix */
13355c5985e7SKris Buschelman     ierr = MatCreate(A->comm,nrows,ncols,nrows,ncols,&newmat);CHKERRQ(ierr);
13365c5985e7SKris Buschelman     ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr);
13375c5985e7SKris Buschelman     ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr);
1338182d2002SSatish Balay   }
1339182d2002SSatish Balay 
1340182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1341182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1342182d2002SSatish Balay 
1343182d2002SSatish Balay   for (i=0; i<ncols; i++) {
1344182d2002SSatish Balay     av = v + m*icol[i];
1345182d2002SSatish Balay     for (j=0; j<nrows; j++) {
1346182d2002SSatish Balay       *bv++ = av[irow[j]];
13470754003eSLois Curfman McInnes     }
13480754003eSLois Curfman McInnes   }
1349182d2002SSatish Balay 
1350182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
13516d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13526d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13530754003eSLois Curfman McInnes 
13540754003eSLois Curfman McInnes   /* Free work space */
135578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
135678b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1357182d2002SSatish Balay   *B = newmat;
13583a40ed3dSBarry Smith   PetscFunctionReturn(0);
13590754003eSLois Curfman McInnes }
13600754003eSLois Curfman McInnes 
13614a2ae208SSatish Balay #undef __FUNCT__
13624a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense"
1363268466fbSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1364905e6a2fSBarry Smith {
1365905e6a2fSBarry Smith   int ierr,i;
1366905e6a2fSBarry Smith 
13673a40ed3dSBarry Smith   PetscFunctionBegin;
1368905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1369b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
1370905e6a2fSBarry Smith   }
1371905e6a2fSBarry Smith 
1372905e6a2fSBarry Smith   for (i=0; i<n; i++) {
13736a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1374905e6a2fSBarry Smith   }
13753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1376905e6a2fSBarry Smith }
1377905e6a2fSBarry Smith 
13784a2ae208SSatish Balay #undef __FUNCT__
13794a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense"
1380cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
13814b0e389bSBarry Smith {
13824b0e389bSBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data;
1383a5ce6ee0Svictorle   int          lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr;
13843a40ed3dSBarry Smith 
13853a40ed3dSBarry Smith   PetscFunctionBegin;
138633f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
138733f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1388cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
13893a40ed3dSBarry Smith     PetscFunctionReturn(0);
13903a40ed3dSBarry Smith   }
13910dbb7854Svictorle   if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1392a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
13930dbb7854Svictorle     for (j=0; j<n; j++) {
1394a5ce6ee0Svictorle       ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr);
1395a5ce6ee0Svictorle     }
1396a5ce6ee0Svictorle   } else {
139787828ca2SBarry Smith     ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr);
1398a5ce6ee0Svictorle   }
1399273d9f13SBarry Smith   PetscFunctionReturn(0);
1400273d9f13SBarry Smith }
1401273d9f13SBarry Smith 
14024a2ae208SSatish Balay #undef __FUNCT__
14034a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense"
1404273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A)
1405273d9f13SBarry Smith {
1406273d9f13SBarry Smith   int        ierr;
1407273d9f13SBarry Smith 
1408273d9f13SBarry Smith   PetscFunctionBegin;
1409273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
14103a40ed3dSBarry Smith   PetscFunctionReturn(0);
14114b0e389bSBarry Smith }
14124b0e389bSBarry Smith 
1413289bc588SBarry Smith /* -------------------------------------------------------------------*/
1414a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense,
1415905e6a2fSBarry Smith        MatGetRow_SeqDense,
1416905e6a2fSBarry Smith        MatRestoreRow_SeqDense,
1417905e6a2fSBarry Smith        MatMult_SeqDense,
141897304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense,
14197c922b88SBarry Smith        MatMultTranspose_SeqDense,
14207c922b88SBarry Smith        MatMultTransposeAdd_SeqDense,
1421905e6a2fSBarry Smith        MatSolve_SeqDense,
1422905e6a2fSBarry Smith        MatSolveAdd_SeqDense,
14237c922b88SBarry Smith        MatSolveTranspose_SeqDense,
142497304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense,
1425905e6a2fSBarry Smith        MatLUFactor_SeqDense,
1426905e6a2fSBarry Smith        MatCholeskyFactor_SeqDense,
1427ec8511deSBarry Smith        MatRelax_SeqDense,
1428ec8511deSBarry Smith        MatTranspose_SeqDense,
142997304618SKris Buschelman /*15*/ MatGetInfo_SeqDense,
1430905e6a2fSBarry Smith        MatEqual_SeqDense,
1431905e6a2fSBarry Smith        MatGetDiagonal_SeqDense,
1432905e6a2fSBarry Smith        MatDiagonalScale_SeqDense,
1433905e6a2fSBarry Smith        MatNorm_SeqDense,
143497304618SKris Buschelman /*20*/ 0,
1435905e6a2fSBarry Smith        0,
1436905e6a2fSBarry Smith        0,
1437905e6a2fSBarry Smith        MatSetOption_SeqDense,
1438905e6a2fSBarry Smith        MatZeroEntries_SeqDense,
143997304618SKris Buschelman /*25*/ MatZeroRows_SeqDense,
1440905e6a2fSBarry Smith        MatLUFactorSymbolic_SeqDense,
1441905e6a2fSBarry Smith        MatLUFactorNumeric_SeqDense,
1442905e6a2fSBarry Smith        MatCholeskyFactorSymbolic_SeqDense,
1443905e6a2fSBarry Smith        MatCholeskyFactorNumeric_SeqDense,
144497304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense,
1445273d9f13SBarry Smith        0,
1446905e6a2fSBarry Smith        0,
1447905e6a2fSBarry Smith        MatGetArray_SeqDense,
1448905e6a2fSBarry Smith        MatRestoreArray_SeqDense,
144997304618SKris Buschelman /*35*/ MatDuplicate_SeqDense,
1450a5ae1ecdSBarry Smith        0,
1451a5ae1ecdSBarry Smith        0,
1452a5ae1ecdSBarry Smith        0,
1453a5ae1ecdSBarry Smith        0,
145497304618SKris Buschelman /*40*/ MatAXPY_SeqDense,
1455a5ae1ecdSBarry Smith        MatGetSubMatrices_SeqDense,
1456a5ae1ecdSBarry Smith        0,
14574b0e389bSBarry Smith        MatGetValues_SeqDense,
1458a5ae1ecdSBarry Smith        MatCopy_SeqDense,
145997304618SKris Buschelman /*45*/ 0,
1460a5ae1ecdSBarry Smith        MatScale_SeqDense,
1461a5ae1ecdSBarry Smith        0,
1462a5ae1ecdSBarry Smith        0,
1463a5ae1ecdSBarry Smith        0,
146497304618SKris Buschelman /*50*/ MatGetBlockSize_SeqDense,
1465a5ae1ecdSBarry Smith        0,
1466a5ae1ecdSBarry Smith        0,
1467a5ae1ecdSBarry Smith        0,
1468a5ae1ecdSBarry Smith        0,
146997304618SKris Buschelman /*55*/ 0,
1470a5ae1ecdSBarry Smith        0,
1471a5ae1ecdSBarry Smith        0,
1472a5ae1ecdSBarry Smith        0,
1473a5ae1ecdSBarry Smith        0,
147497304618SKris Buschelman /*60*/ 0,
1475e03a110bSBarry Smith        MatDestroy_SeqDense,
1476e03a110bSBarry Smith        MatView_SeqDense,
147797304618SKris Buschelman        MatGetPetscMaps_Petsc,
147897304618SKris Buschelman        0,
147997304618SKris Buschelman /*65*/ 0,
148097304618SKris Buschelman        0,
148197304618SKris Buschelman        0,
148297304618SKris Buschelman        0,
148397304618SKris Buschelman        0,
148497304618SKris Buschelman /*70*/ 0,
148597304618SKris Buschelman        0,
148697304618SKris Buschelman        0,
148797304618SKris Buschelman        0,
148897304618SKris Buschelman        0,
148997304618SKris Buschelman /*75*/ 0,
149097304618SKris Buschelman        0,
149197304618SKris Buschelman        0,
149297304618SKris Buschelman        0,
149397304618SKris Buschelman        0,
149497304618SKris Buschelman /*80*/ 0,
149597304618SKris Buschelman        0,
149697304618SKris Buschelman        0,
149797304618SKris Buschelman        0,
149897304618SKris Buschelman /*85*/ MatLoad_SeqDense};
149990ace30eSBarry Smith 
15004a2ae208SSatish Balay #undef __FUNCT__
15014a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense"
15024b828684SBarry Smith /*@C
1503fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
1504d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
1505d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
1506289bc588SBarry Smith 
1507db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1508db81eaa0SLois Curfman McInnes 
150920563c6bSBarry Smith    Input Parameters:
1510db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
15110c775827SLois Curfman McInnes .  m - number of rows
151218f449edSLois Curfman McInnes .  n - number of columns
1513db81eaa0SLois Curfman McInnes -  data - optional location of matrix data.  Set data=PETSC_NULL for PETSc
1514dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
151520563c6bSBarry Smith 
151620563c6bSBarry Smith    Output Parameter:
151744cd7ae7SLois Curfman McInnes .  A - the matrix
151820563c6bSBarry Smith 
1519b259b22eSLois Curfman McInnes    Notes:
152018f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
152118f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
1522b4fd4287SBarry Smith    set data=PETSC_NULL.
152318f449edSLois Curfman McInnes 
1524027ccd11SLois Curfman McInnes    Level: intermediate
1525027ccd11SLois Curfman McInnes 
1526dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS
1527d65003e9SLois Curfman McInnes 
1528db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
152920563c6bSBarry Smith @*/
153087828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A)
1531289bc588SBarry Smith {
1532273d9f13SBarry Smith   int ierr;
15333b2fbd54SBarry Smith 
15343a40ed3dSBarry Smith   PetscFunctionBegin;
1535273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1536273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
1537273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
1538273d9f13SBarry Smith   PetscFunctionReturn(0);
1539273d9f13SBarry Smith }
1540273d9f13SBarry Smith 
15414a2ae208SSatish Balay #undef __FUNCT__
15424a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation"
1543273d9f13SBarry Smith /*@C
1544273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
1545273d9f13SBarry Smith 
1546273d9f13SBarry Smith    Collective on MPI_Comm
1547273d9f13SBarry Smith 
1548273d9f13SBarry Smith    Input Parameters:
1549273d9f13SBarry Smith +  A - the matrix
1550273d9f13SBarry Smith -  data - the array (or PETSC_NULL)
1551273d9f13SBarry Smith 
1552273d9f13SBarry Smith    Notes:
1553273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
1554273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
1555273d9f13SBarry Smith    set data=PETSC_NULL.
1556273d9f13SBarry Smith 
1557273d9f13SBarry Smith    Level: intermediate
1558273d9f13SBarry Smith 
1559273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS
1560273d9f13SBarry Smith 
1561273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues()
1562273d9f13SBarry Smith @*/
1563ca01db9bSBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
1564273d9f13SBarry Smith {
1565ca01db9bSBarry Smith   int ierr,(*f)(Mat,PetscScalar[]);
1566a23d5eceSKris Buschelman 
1567a23d5eceSKris Buschelman   PetscFunctionBegin;
1568a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1569a23d5eceSKris Buschelman   if (f) {
1570a23d5eceSKris Buschelman     ierr = (*f)(B,data);CHKERRQ(ierr);
1571a23d5eceSKris Buschelman   }
1572a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1573a23d5eceSKris Buschelman }
1574a23d5eceSKris Buschelman 
1575a23d5eceSKris Buschelman EXTERN_C_BEGIN
1576a23d5eceSKris Buschelman #undef __FUNCT__
1577a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense"
1578a23d5eceSKris Buschelman int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
1579a23d5eceSKris Buschelman {
1580273d9f13SBarry Smith   Mat_SeqDense *b;
1581273d9f13SBarry Smith   int          ierr;
1582273d9f13SBarry Smith 
1583273d9f13SBarry Smith   PetscFunctionBegin;
1584273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1585273d9f13SBarry Smith   b               = (Mat_SeqDense*)B->data;
1586273d9f13SBarry Smith   if (!data) {
158787828ca2SBarry Smith     ierr          = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr);
158887828ca2SBarry Smith     ierr          = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr);
1589273d9f13SBarry Smith     b->user_alloc = PETSC_FALSE;
159087828ca2SBarry Smith     PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));
1591273d9f13SBarry Smith   } else { /* user-allocated storage */
1592273d9f13SBarry Smith     b->v          = data;
1593273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
1594273d9f13SBarry Smith   }
1595273d9f13SBarry Smith   PetscFunctionReturn(0);
1596273d9f13SBarry Smith }
1597a23d5eceSKris Buschelman EXTERN_C_END
1598273d9f13SBarry Smith 
15991b807ce4Svictorle #undef __FUNCT__
16001b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA"
16011b807ce4Svictorle /*@C
16021b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
16031b807ce4Svictorle 
16041b807ce4Svictorle   Input parameter:
16051b807ce4Svictorle + A - the matrix
16061b807ce4Svictorle - lda - the leading dimension
16071b807ce4Svictorle 
16081b807ce4Svictorle   Notes:
16091b807ce4Svictorle   This routine is to be used in conjunction with MatSeqDenseSetPreallocation;
16101b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
16111b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
16121b807ce4Svictorle 
16131b807ce4Svictorle   Level: intermediate
16141b807ce4Svictorle 
16151b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS
16161b807ce4Svictorle 
16171b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation()
16181b807ce4Svictorle @*/
16191b807ce4Svictorle int MatSeqDenseSetLDA(Mat B,int lda)
16201b807ce4Svictorle {
16211b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
16221b807ce4Svictorle   PetscFunctionBegin;
1623*634064b4SBarry Smith   if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %d must be at least matrix dimension %d",lda,B->m);
16241b807ce4Svictorle   b->lda = lda;
16251b807ce4Svictorle   PetscFunctionReturn(0);
16261b807ce4Svictorle }
16271b807ce4Svictorle 
16280bad9183SKris Buschelman /*MC
1629fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
16300bad9183SKris Buschelman 
16310bad9183SKris Buschelman    Options Database Keys:
16320bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
16330bad9183SKris Buschelman 
16340bad9183SKris Buschelman   Level: beginner
16350bad9183SKris Buschelman 
16360bad9183SKris Buschelman .seealso: MatCreateSeqDense
16370bad9183SKris Buschelman M*/
16380bad9183SKris Buschelman 
1639273d9f13SBarry Smith EXTERN_C_BEGIN
16404a2ae208SSatish Balay #undef __FUNCT__
16414a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense"
1642273d9f13SBarry Smith int MatCreate_SeqDense(Mat B)
1643273d9f13SBarry Smith {
1644273d9f13SBarry Smith   Mat_SeqDense *b;
1645273d9f13SBarry Smith   int          ierr,size;
1646273d9f13SBarry Smith 
1647273d9f13SBarry Smith   PetscFunctionBegin;
1648273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
164929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
165055659b69SBarry Smith 
1651273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1652273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1653273d9f13SBarry Smith 
1654b0a32e0cSBarry Smith   ierr            = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr);
1655549d3d68SSatish Balay   ierr            = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr);
1656549d3d68SSatish Balay   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
165744cd7ae7SLois Curfman McInnes   B->factor       = 0;
165890f02eecSBarry Smith   B->mapping      = 0;
1659b0a32e0cSBarry Smith   PetscLogObjectMemory(B,sizeof(struct _p_Mat));
166044cd7ae7SLois Curfman McInnes   B->data         = (void*)b;
166118f449edSLois Curfman McInnes 
16628a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
16638a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1664a5ae1ecdSBarry Smith 
166544cd7ae7SLois Curfman McInnes   b->pivots       = 0;
1666273d9f13SBarry Smith   b->roworiented  = PETSC_TRUE;
1667273d9f13SBarry Smith   b->v            = 0;
16681b807ce4Svictorle   b->lda          = B->m;
16694e220ebcSLois Curfman McInnes 
1670a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C",
1671a23d5eceSKris Buschelman                                     "MatSeqDenseSetPreallocation_SeqDense",
1672a23d5eceSKris Buschelman                                      MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
16733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1674289bc588SBarry Smith }
1675273d9f13SBarry Smith EXTERN_C_END
1676