xref: /petsc/src/mat/impls/dense/seq/dense.c (revision 071fcb05f2a6aea8aef2c2090580530b9e9df77e)
1be1d678aSKris Buschelman 
267e560aaSBarry Smith /*
367e560aaSBarry Smith      Defines the basic matrix operations for sequential dense.
467e560aaSBarry Smith */
5289bc588SBarry Smith 
6dec5eb66SMatthew G Knepley #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <petscblaslapack.h>
8289bc588SBarry Smith 
96a63e612SBarry Smith #include <../src/mat/impls/aij/seq/aij.h>
10b2573a8aSBarry Smith 
118c178816SStefano Zampini static PetscErrorCode MatSeqDenseSymmetrize_Private(Mat A, PetscBool hermitian)
128c178816SStefano Zampini {
138c178816SStefano Zampini   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
148c178816SStefano Zampini   PetscInt      j, k, n = A->rmap->n;
158c178816SStefano Zampini 
168c178816SStefano Zampini   PetscFunctionBegin;
178c178816SStefano Zampini   if (A->rmap->n != A->cmap->n) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot symmetrize a rectangular matrix");
188c178816SStefano Zampini   if (!hermitian) {
198c178816SStefano Zampini     for (k=0;k<n;k++) {
208c178816SStefano Zampini       for (j=k;j<n;j++) {
218c178816SStefano Zampini         mat->v[j*mat->lda + k] = mat->v[k*mat->lda + j];
228c178816SStefano Zampini       }
238c178816SStefano Zampini     }
248c178816SStefano Zampini   } else {
258c178816SStefano Zampini     for (k=0;k<n;k++) {
268c178816SStefano Zampini       for (j=k;j<n;j++) {
278c178816SStefano Zampini         mat->v[j*mat->lda + k] = PetscConj(mat->v[k*mat->lda + j]);
288c178816SStefano Zampini       }
298c178816SStefano Zampini     }
308c178816SStefano Zampini   }
318c178816SStefano Zampini   PetscFunctionReturn(0);
328c178816SStefano Zampini }
338c178816SStefano Zampini 
3405709791SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat A)
358c178816SStefano Zampini {
368c178816SStefano Zampini #if defined(PETSC_MISSING_LAPACK_POTRF)
378c178816SStefano Zampini   PetscFunctionBegin;
388c178816SStefano Zampini   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
398c178816SStefano Zampini #else
408c178816SStefano Zampini   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
418c178816SStefano Zampini   PetscErrorCode ierr;
428c178816SStefano Zampini   PetscBLASInt   info,n;
438c178816SStefano Zampini 
448c178816SStefano Zampini   PetscFunctionBegin;
458c178816SStefano Zampini   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
468c178816SStefano Zampini   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
478c178816SStefano Zampini   if (A->factortype == MAT_FACTOR_LU) {
488c178816SStefano Zampini     if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
498c178816SStefano Zampini     if (!mat->fwork) {
508c178816SStefano Zampini       mat->lfwork = n;
518c178816SStefano Zampini       ierr = PetscMalloc1(mat->lfwork,&mat->fwork);CHKERRQ(ierr);
528c178816SStefano Zampini       ierr = PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));CHKERRQ(ierr);
538c178816SStefano Zampini     }
5400121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
558c178816SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
5600121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
578c178816SStefano Zampini     ierr = PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); /* TODO CHECK FLOPS */
588c178816SStefano Zampini   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
598c178816SStefano Zampini     if (A->spd) {
6000121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
618c178816SStefano Zampini       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info));
6200121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
638c178816SStefano Zampini       ierr = MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);CHKERRQ(ierr);
648c178816SStefano Zampini #if defined(PETSC_USE_COMPLEX)
658c178816SStefano Zampini     } else if (A->hermitian) {
668c178816SStefano Zampini       if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
678c178816SStefano Zampini       if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
6800121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
698c178816SStefano Zampini       PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
7000121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
718c178816SStefano Zampini       ierr = MatSeqDenseSymmetrize_Private(A,PETSC_TRUE);CHKERRQ(ierr);
728c178816SStefano Zampini #endif
738c178816SStefano Zampini     } else { /* symmetric case */
748c178816SStefano Zampini       if (!mat->pivots) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present");
758c178816SStefano Zampini       if (!mat->fwork) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present");
7600121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
778c178816SStefano Zampini       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info));
7800121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
798c178816SStefano Zampini       ierr = MatSeqDenseSymmetrize_Private(A,PETSC_FALSE);CHKERRQ(ierr);
808c178816SStefano Zampini     }
818c178816SStefano Zampini     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %D",(PetscInt)info-1);
828c178816SStefano Zampini     ierr = PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); /* TODO CHECK FLOPS */
838c178816SStefano Zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
848c178816SStefano Zampini #endif
858c178816SStefano Zampini 
868c178816SStefano Zampini   A->ops->solve             = NULL;
878c178816SStefano Zampini   A->ops->matsolve          = NULL;
888c178816SStefano Zampini   A->ops->solvetranspose    = NULL;
898c178816SStefano Zampini   A->ops->matsolvetranspose = NULL;
908c178816SStefano Zampini   A->ops->solveadd          = NULL;
918c178816SStefano Zampini   A->ops->solvetransposeadd = NULL;
928c178816SStefano Zampini   A->factortype             = MAT_FACTOR_NONE;
938c178816SStefano Zampini   ierr                      = PetscFree(A->solvertype);CHKERRQ(ierr);
948c178816SStefano Zampini   PetscFunctionReturn(0);
958c178816SStefano Zampini }
968c178816SStefano Zampini 
973f49a652SStefano Zampini PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
983f49a652SStefano Zampini {
993f49a652SStefano Zampini   PetscErrorCode    ierr;
1003f49a652SStefano Zampini   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1013f49a652SStefano Zampini   PetscInt          m  = l->lda, n = A->cmap->n,r = A->rmap->n, i,j;
1023f49a652SStefano Zampini   PetscScalar       *slot,*bb;
1033f49a652SStefano Zampini   const PetscScalar *xx;
1043f49a652SStefano Zampini 
1053f49a652SStefano Zampini   PetscFunctionBegin;
1063f49a652SStefano Zampini #if defined(PETSC_USE_DEBUG)
1073f49a652SStefano Zampini   for (i=0; i<N; i++) {
1083f49a652SStefano Zampini     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1093f49a652SStefano Zampini     if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
1103f49a652SStefano Zampini     if (rows[i] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Col %D requested to be zeroed greater than or equal number of cols %D",rows[i],A->cmap->n);
1113f49a652SStefano Zampini   }
1123f49a652SStefano Zampini #endif
1133f49a652SStefano Zampini 
1143f49a652SStefano Zampini   /* fix right hand side if needed */
1153f49a652SStefano Zampini   if (x && b) {
1166c4d906cSStefano Zampini     Vec xt;
1176c4d906cSStefano Zampini 
1186c4d906cSStefano Zampini     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1196c4d906cSStefano Zampini     ierr = VecDuplicate(x,&xt);CHKERRQ(ierr);
1206c4d906cSStefano Zampini     ierr = VecCopy(x,xt);CHKERRQ(ierr);
1216c4d906cSStefano Zampini     ierr = VecScale(xt,-1.0);CHKERRQ(ierr);
1226c4d906cSStefano Zampini     ierr = MatMultAdd(A,xt,b,b);CHKERRQ(ierr);
1236c4d906cSStefano Zampini     ierr = VecDestroy(&xt);CHKERRQ(ierr);
1243f49a652SStefano Zampini     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1253f49a652SStefano Zampini     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
1263f49a652SStefano Zampini     for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
1273f49a652SStefano Zampini     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1283f49a652SStefano Zampini     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
1293f49a652SStefano Zampini   }
1303f49a652SStefano Zampini 
1313f49a652SStefano Zampini   for (i=0; i<N; i++) {
1323f49a652SStefano Zampini     slot = l->v + rows[i]*m;
133580bdb30SBarry Smith     ierr = PetscArrayzero(slot,r);CHKERRQ(ierr);
1343f49a652SStefano Zampini   }
1353f49a652SStefano Zampini   for (i=0; i<N; i++) {
1363f49a652SStefano Zampini     slot = l->v + rows[i];
1373f49a652SStefano Zampini     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
1383f49a652SStefano Zampini   }
1393f49a652SStefano Zampini   if (diag != 0.0) {
1403f49a652SStefano Zampini     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
1413f49a652SStefano Zampini     for (i=0; i<N; i++) {
1423f49a652SStefano Zampini       slot  = l->v + (m+1)*rows[i];
1433f49a652SStefano Zampini       *slot = diag;
1443f49a652SStefano Zampini     }
1453f49a652SStefano Zampini   }
1463f49a652SStefano Zampini   PetscFunctionReturn(0);
1473f49a652SStefano Zampini }
1483f49a652SStefano Zampini 
149abc3b08eSStefano Zampini PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C)
150abc3b08eSStefano Zampini {
151abc3b08eSStefano Zampini   Mat_SeqDense   *c = (Mat_SeqDense*)(C->data);
152abc3b08eSStefano Zampini   PetscErrorCode ierr;
153abc3b08eSStefano Zampini 
154abc3b08eSStefano Zampini   PetscFunctionBegin;
155abc3b08eSStefano Zampini   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,P,c->ptapwork);CHKERRQ(ierr);
156abc3b08eSStefano Zampini   ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(P,c->ptapwork,C);CHKERRQ(ierr);
157abc3b08eSStefano Zampini   PetscFunctionReturn(0);
158abc3b08eSStefano Zampini }
159abc3b08eSStefano Zampini 
160abc3b08eSStefano Zampini PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat *C)
161abc3b08eSStefano Zampini {
162abc3b08eSStefano Zampini   Mat_SeqDense   *c;
163abc3b08eSStefano Zampini   PetscErrorCode ierr;
164abc3b08eSStefano Zampini 
165abc3b08eSStefano Zampini   PetscFunctionBegin;
166abc3b08eSStefano Zampini   ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)A),P->cmap->N,P->cmap->N,NULL,C);CHKERRQ(ierr);
167abc3b08eSStefano Zampini   c = (Mat_SeqDense*)((*C)->data);
168abc3b08eSStefano Zampini   ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)A),A->rmap->N,P->cmap->N,NULL,&c->ptapwork);CHKERRQ(ierr);
169abc3b08eSStefano Zampini   PetscFunctionReturn(0);
170abc3b08eSStefano Zampini }
171abc3b08eSStefano Zampini 
172150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_SeqDense_SeqDense(Mat A,Mat P,MatReuse reuse,PetscReal fill,Mat *C)
173abc3b08eSStefano Zampini {
174abc3b08eSStefano Zampini   PetscErrorCode ierr;
175abc3b08eSStefano Zampini 
176abc3b08eSStefano Zampini   PetscFunctionBegin;
177abc3b08eSStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
178abc3b08eSStefano Zampini     ierr = MatPtAPSymbolic_SeqDense_SeqDense(A,P,fill,C);CHKERRQ(ierr);
179abc3b08eSStefano Zampini   }
180abc3b08eSStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
181abc3b08eSStefano Zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
182abc3b08eSStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
183abc3b08eSStefano Zampini   PetscFunctionReturn(0);
184abc3b08eSStefano Zampini }
185abc3b08eSStefano Zampini 
186cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
187b49cda9fSStefano Zampini {
188a13144ffSStefano Zampini   Mat            B = NULL;
189b49cda9fSStefano Zampini   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
190b49cda9fSStefano Zampini   Mat_SeqDense   *b;
191b49cda9fSStefano Zampini   PetscErrorCode ierr;
192b49cda9fSStefano Zampini   PetscInt       *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i;
193b49cda9fSStefano Zampini   MatScalar      *av=a->a;
194a13144ffSStefano Zampini   PetscBool      isseqdense;
195b49cda9fSStefano Zampini 
196b49cda9fSStefano Zampini   PetscFunctionBegin;
197a13144ffSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
198a13144ffSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense);CHKERRQ(ierr);
199a13144ffSStefano Zampini     if (!isseqdense) SETERRQ1(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type);
200a13144ffSStefano Zampini   }
201a13144ffSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
202b49cda9fSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
203b49cda9fSStefano Zampini     ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
204b49cda9fSStefano Zampini     ierr = MatSetType(B,MATSEQDENSE);CHKERRQ(ierr);
205b49cda9fSStefano Zampini     ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
206b49cda9fSStefano Zampini     b    = (Mat_SeqDense*)(B->data);
207a13144ffSStefano Zampini   } else {
208a13144ffSStefano Zampini     b    = (Mat_SeqDense*)((*newmat)->data);
209580bdb30SBarry Smith     ierr = PetscArrayzero(b->v,m*n);CHKERRQ(ierr);
210a13144ffSStefano Zampini   }
211b49cda9fSStefano Zampini   for (i=0; i<m; i++) {
212b49cda9fSStefano Zampini     PetscInt j;
213b49cda9fSStefano Zampini     for (j=0;j<ai[1]-ai[0];j++) {
214b49cda9fSStefano Zampini       b->v[*aj*m+i] = *av;
215b49cda9fSStefano Zampini       aj++;
216b49cda9fSStefano Zampini       av++;
217b49cda9fSStefano Zampini     }
218b49cda9fSStefano Zampini     ai++;
219b49cda9fSStefano Zampini   }
220b49cda9fSStefano Zampini 
221511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
222a13144ffSStefano Zampini     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223a13144ffSStefano Zampini     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
22428be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
225b49cda9fSStefano Zampini   } else {
226a13144ffSStefano Zampini     if (B) *newmat = B;
227a13144ffSStefano Zampini     ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
228a13144ffSStefano Zampini     ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
229b49cda9fSStefano Zampini   }
230b49cda9fSStefano Zampini   PetscFunctionReturn(0);
231b49cda9fSStefano Zampini }
232b49cda9fSStefano Zampini 
233cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
2346a63e612SBarry Smith {
2356a63e612SBarry Smith   Mat            B;
2366a63e612SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2376a63e612SBarry Smith   PetscErrorCode ierr;
2389399e1b8SMatthew G. Knepley   PetscInt       i, j;
2399399e1b8SMatthew G. Knepley   PetscInt       *rows, *nnz;
2409399e1b8SMatthew G. Knepley   MatScalar      *aa = a->v, *vals;
2416a63e612SBarry Smith 
2426a63e612SBarry Smith   PetscFunctionBegin;
243ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
2446a63e612SBarry Smith   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2456a63e612SBarry Smith   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
2469399e1b8SMatthew G. Knepley   ierr = PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);CHKERRQ(ierr);
2479399e1b8SMatthew G. Knepley   for (j=0; j<A->cmap->n; j++) {
2489399e1b8SMatthew G. Knepley     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i];
2496a63e612SBarry Smith     aa += a->lda;
2506a63e612SBarry Smith   }
2519399e1b8SMatthew G. Knepley   ierr = MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);CHKERRQ(ierr);
2529399e1b8SMatthew G. Knepley   aa = a->v;
2539399e1b8SMatthew G. Knepley   for (j=0; j<A->cmap->n; j++) {
2549399e1b8SMatthew G. Knepley     PetscInt numRows = 0;
2559399e1b8SMatthew G. Knepley     for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];}
2569399e1b8SMatthew G. Knepley     ierr = MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);CHKERRQ(ierr);
2579399e1b8SMatthew G. Knepley     aa  += a->lda;
2589399e1b8SMatthew G. Knepley   }
2599399e1b8SMatthew G. Knepley   ierr = PetscFree3(rows,nnz,vals);CHKERRQ(ierr);
2606a63e612SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2616a63e612SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2626a63e612SBarry Smith 
263511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
26428be2f97SBarry Smith     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
2656a63e612SBarry Smith   } else {
2666a63e612SBarry Smith     *newmat = B;
2676a63e612SBarry Smith   }
2686a63e612SBarry Smith   PetscFunctionReturn(0);
2696a63e612SBarry Smith }
2706a63e612SBarry Smith 
271e0877f53SBarry Smith static PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str)
2721987afe7SBarry Smith {
2731987afe7SBarry Smith   Mat_SeqDense   *x     = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data;
274f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
27513f74950SBarry Smith   PetscInt       j;
2760805154bSBarry Smith   PetscBLASInt   N,m,ldax,lday,one = 1;
277efee365bSSatish Balay   PetscErrorCode ierr;
2783a40ed3dSBarry Smith 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
280c5df96a5SBarry Smith   ierr = PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);CHKERRQ(ierr);
281c5df96a5SBarry Smith   ierr = PetscBLASIntCast(X->rmap->n,&m);CHKERRQ(ierr);
282c5df96a5SBarry Smith   ierr = PetscBLASIntCast(x->lda,&ldax);CHKERRQ(ierr);
283c5df96a5SBarry Smith   ierr = PetscBLASIntCast(y->lda,&lday);CHKERRQ(ierr);
284a5ce6ee0Svictorle   if (ldax>m || lday>m) {
285d0f46423SBarry Smith     for (j=0; j<X->cmap->n; j++) {
2868b83055fSJed Brown       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one));
287a5ce6ee0Svictorle     }
288a5ce6ee0Svictorle   } else {
2898b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one));
290a5ce6ee0Svictorle   }
291a3fa217bSJose E. Roman   ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
2920450473dSBarry Smith   ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr);
2933a40ed3dSBarry Smith   PetscFunctionReturn(0);
2941987afe7SBarry Smith }
2951987afe7SBarry Smith 
296e0877f53SBarry Smith static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info)
297289bc588SBarry Smith {
298d0f46423SBarry Smith   PetscInt N = A->rmap->n*A->cmap->n;
2993a40ed3dSBarry Smith 
3003a40ed3dSBarry Smith   PetscFunctionBegin;
3014e220ebcSLois Curfman McInnes   info->block_size        = 1.0;
3024e220ebcSLois Curfman McInnes   info->nz_allocated      = (double)N;
3036de62eeeSBarry Smith   info->nz_used           = (double)N;
3046de62eeeSBarry Smith   info->nz_unneeded       = (double)0;
3054e220ebcSLois Curfman McInnes   info->assemblies        = (double)A->num_ass;
3064e220ebcSLois Curfman McInnes   info->mallocs           = 0;
3077adad957SLisandro Dalcin   info->memory            = ((PetscObject)A)->mem;
3084e220ebcSLois Curfman McInnes   info->fill_ratio_given  = 0;
3094e220ebcSLois Curfman McInnes   info->fill_ratio_needed = 0;
3104e220ebcSLois Curfman McInnes   info->factor_mallocs    = 0;
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
312289bc588SBarry Smith }
313289bc588SBarry Smith 
314e0877f53SBarry Smith static PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha)
31580cd9d93SLois Curfman McInnes {
316273d9f13SBarry Smith   Mat_SeqDense   *a     = (Mat_SeqDense*)A->data;
317f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
318efee365bSSatish Balay   PetscErrorCode ierr;
319c5df96a5SBarry Smith   PetscBLASInt   one = 1,j,nz,lda;
32080cd9d93SLois Curfman McInnes 
3213a40ed3dSBarry Smith   PetscFunctionBegin;
322c5df96a5SBarry Smith   ierr = PetscBLASIntCast(a->lda,&lda);CHKERRQ(ierr);
323d0f46423SBarry Smith   if (lda>A->rmap->n) {
324c5df96a5SBarry Smith     ierr = PetscBLASIntCast(A->rmap->n,&nz);CHKERRQ(ierr);
325d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
3268b83055fSJed Brown       PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v+j*lda,&one));
327a5ce6ee0Svictorle     }
328a5ce6ee0Svictorle   } else {
329c5df96a5SBarry Smith     ierr = PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);CHKERRQ(ierr);
3308b83055fSJed Brown     PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v,&one));
331a5ce6ee0Svictorle   }
332efee365bSSatish Balay   ierr = PetscLogFlops(nz);CHKERRQ(ierr);
3333a40ed3dSBarry Smith   PetscFunctionReturn(0);
33480cd9d93SLois Curfman McInnes }
33580cd9d93SLois Curfman McInnes 
336e0877f53SBarry Smith static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool  *fl)
3371cbb95d3SBarry Smith {
3381cbb95d3SBarry Smith   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
339d0f46423SBarry Smith   PetscInt     i,j,m = A->rmap->n,N;
3401cbb95d3SBarry Smith   PetscScalar  *v = a->v;
3411cbb95d3SBarry Smith 
3421cbb95d3SBarry Smith   PetscFunctionBegin;
3431cbb95d3SBarry Smith   *fl = PETSC_FALSE;
344d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0);
3451cbb95d3SBarry Smith   N = a->lda;
3461cbb95d3SBarry Smith 
3471cbb95d3SBarry Smith   for (i=0; i<m; i++) {
3481cbb95d3SBarry Smith     for (j=i+1; j<m; j++) {
3491cbb95d3SBarry Smith       if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0);
3501cbb95d3SBarry Smith     }
3511cbb95d3SBarry Smith   }
3521cbb95d3SBarry Smith   *fl = PETSC_TRUE;
3531cbb95d3SBarry Smith   PetscFunctionReturn(0);
3541cbb95d3SBarry Smith }
3551cbb95d3SBarry Smith 
356e0877f53SBarry Smith static PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues)
357b24902e0SBarry Smith {
358b24902e0SBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data,*l;
359b24902e0SBarry Smith   PetscErrorCode ierr;
360b24902e0SBarry Smith   PetscInt       lda = (PetscInt)mat->lda,j,m;
361b24902e0SBarry Smith 
362b24902e0SBarry Smith   PetscFunctionBegin;
363aa5ea44dSBarry Smith   ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr);
364aa5ea44dSBarry Smith   ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr);
3650298fd71SBarry Smith   ierr = MatSeqDenseSetPreallocation(newi,NULL);CHKERRQ(ierr);
366b24902e0SBarry Smith   if (cpvalues == MAT_COPY_VALUES) {
367b24902e0SBarry Smith     l = (Mat_SeqDense*)newi->data;
368d0f46423SBarry Smith     if (lda>A->rmap->n) {
369d0f46423SBarry Smith       m = A->rmap->n;
370d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
371580bdb30SBarry Smith         ierr = PetscArraycpy(l->v+j*m,mat->v+j*lda,m);CHKERRQ(ierr);
372b24902e0SBarry Smith       }
373b24902e0SBarry Smith     } else {
374580bdb30SBarry Smith       ierr = PetscArraycpy(l->v,mat->v,A->rmap->n*A->cmap->n);CHKERRQ(ierr);
375b24902e0SBarry Smith     }
376b24902e0SBarry Smith   }
377b24902e0SBarry Smith   newi->assembled = PETSC_TRUE;
378b24902e0SBarry Smith   PetscFunctionReturn(0);
379b24902e0SBarry Smith }
380b24902e0SBarry Smith 
381e0877f53SBarry Smith static PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat)
38202cad45dSBarry Smith {
3836849ba73SBarry Smith   PetscErrorCode ierr;
38402cad45dSBarry Smith 
3853a40ed3dSBarry Smith   PetscFunctionBegin;
386ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),newmat);CHKERRQ(ierr);
387d0f46423SBarry Smith   ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
3885c9eb25fSBarry Smith   ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
389719d5645SBarry Smith   ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr);
390b24902e0SBarry Smith   PetscFunctionReturn(0);
391b24902e0SBarry Smith }
392b24902e0SBarry Smith 
3936ee01492SSatish Balay 
394e0877f53SBarry Smith static PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
395719d5645SBarry Smith 
396e0877f53SBarry Smith static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
397289bc588SBarry Smith {
3984482741eSBarry Smith   MatFactorInfo  info;
399a093e273SMatthew Knepley   PetscErrorCode ierr;
4003a40ed3dSBarry Smith 
4013a40ed3dSBarry Smith   PetscFunctionBegin;
402c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
403719d5645SBarry Smith   ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr);
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
405289bc588SBarry Smith }
4066ee01492SSatish Balay 
407e0877f53SBarry Smith static PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy)
408289bc588SBarry Smith {
409c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
4106849ba73SBarry Smith   PetscErrorCode    ierr;
411f1ceaac6SMatthew G. Knepley   const PetscScalar *x;
412f1ceaac6SMatthew G. Knepley   PetscScalar       *y;
413c5df96a5SBarry Smith   PetscBLASInt      one = 1,info,m;
41467e560aaSBarry Smith 
4153a40ed3dSBarry Smith   PetscFunctionBegin;
416c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
417f1ceaac6SMatthew G. Knepley   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
4181ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
419580bdb30SBarry Smith   ierr = PetscArraycpy(y,x,A->rmap->n);CHKERRQ(ierr);
420d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
421ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
422e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
423ae7cfcebSSatish Balay #else
42400121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4258b83055fSJed Brown     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
42600121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
427e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
428ae7cfcebSSatish Balay #endif
429d5f3da31SBarry Smith   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
430ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
431e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
432ae7cfcebSSatish Balay #else
433a49dc2a2SStefano Zampini     if (A->spd) {
43400121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4358b83055fSJed Brown       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
43600121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
437e32f2f54SBarry Smith       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
438a49dc2a2SStefano Zampini #if defined(PETSC_USE_COMPLEX)
439a49dc2a2SStefano Zampini     } else if (A->hermitian) {
44000121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
441a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
44200121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
443a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
444a49dc2a2SStefano Zampini #endif
445a49dc2a2SStefano Zampini     } else { /* symmetric case */
44600121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
447a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
44800121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
449a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
450a49dc2a2SStefano Zampini     }
451ae7cfcebSSatish Balay #endif
4522205254eSKarl Rupp   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
453f1ceaac6SMatthew G. Knepley   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
4541ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
455dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
4563a40ed3dSBarry Smith   PetscFunctionReturn(0);
457289bc588SBarry Smith }
4586ee01492SSatish Balay 
459e0877f53SBarry Smith static PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X)
46085e2c93fSHong Zhang {
46185e2c93fSHong Zhang   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
46285e2c93fSHong Zhang   PetscErrorCode    ierr;
4631683a169SBarry Smith   const PetscScalar *b;
4641683a169SBarry Smith   PetscScalar       *x;
465efb80c78SLisandro Dalcin   PetscInt          n;
466783b601eSJed Brown   PetscBLASInt      nrhs,info,m;
467bda8bf91SBarry Smith   PetscBool         flg;
46885e2c93fSHong Zhang 
46985e2c93fSHong Zhang   PetscFunctionBegin;
470c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
4710298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
472ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
4730298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
474ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
475bda8bf91SBarry Smith 
4760298fd71SBarry Smith   ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr);
477c5df96a5SBarry Smith   ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr);
4781683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr);
4798c778c55SBarry Smith   ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr);
48085e2c93fSHong Zhang 
481580bdb30SBarry Smith   ierr = PetscArraycpy(x,b,m*nrhs);CHKERRQ(ierr);
48285e2c93fSHong Zhang 
48385e2c93fSHong Zhang   if (A->factortype == MAT_FACTOR_LU) {
48485e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS)
48585e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
48685e2c93fSHong Zhang #else
48700121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4888b83055fSJed Brown     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
48900121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
49085e2c93fSHong Zhang     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve");
49185e2c93fSHong Zhang #endif
49285e2c93fSHong Zhang   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
49385e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS)
49485e2c93fSHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
49585e2c93fSHong Zhang #else
496a49dc2a2SStefano Zampini     if (A->spd) {
49700121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
4988b83055fSJed Brown       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info));
49900121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
50085e2c93fSHong Zhang       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
501a49dc2a2SStefano Zampini #if defined(PETSC_USE_COMPLEX)
502a49dc2a2SStefano Zampini     } else if (A->hermitian) {
50300121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
504a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
50500121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
506a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve");
507a49dc2a2SStefano Zampini #endif
508a49dc2a2SStefano Zampini     } else { /* symmetric case */
50900121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
510a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info));
51100121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
512a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
513a49dc2a2SStefano Zampini     }
51485e2c93fSHong Zhang #endif
5152205254eSKarl Rupp   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
51685e2c93fSHong Zhang 
5171683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr);
5188c778c55SBarry Smith   ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr);
51985e2c93fSHong Zhang   ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr);
52085e2c93fSHong Zhang   PetscFunctionReturn(0);
52185e2c93fSHong Zhang }
52285e2c93fSHong Zhang 
52300121966SStefano Zampini static PetscErrorCode MatConjugate_SeqDense(Mat);
52400121966SStefano Zampini 
525e0877f53SBarry Smith static PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy)
526da3a660dSBarry Smith {
527c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
528dfbe8321SBarry Smith   PetscErrorCode    ierr;
529f1ceaac6SMatthew G. Knepley   const PetscScalar *x;
530f1ceaac6SMatthew G. Knepley   PetscScalar       *y;
531c5df96a5SBarry Smith   PetscBLASInt      one = 1,info,m;
53267e560aaSBarry Smith 
5333a40ed3dSBarry Smith   PetscFunctionBegin;
534c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
535f1ceaac6SMatthew G. Knepley   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
5361ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
537580bdb30SBarry Smith   ierr = PetscArraycpy(y,x,A->rmap->n);CHKERRQ(ierr);
5388208b9aeSStefano Zampini   if (A->factortype == MAT_FACTOR_LU) {
539ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS)
540e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable.");
541ae7cfcebSSatish Balay #else
54200121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5438b83055fSJed Brown     PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
54400121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
545e32f2f54SBarry Smith     if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve");
546ae7cfcebSSatish Balay #endif
5478208b9aeSStefano Zampini   } else if (A->factortype == MAT_FACTOR_CHOLESKY) {
548ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS)
549e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable.");
550ae7cfcebSSatish Balay #else
551a49dc2a2SStefano Zampini     if (A->spd) {
55200121966SStefano Zampini #if defined(PETSC_USE_COMPLEX)
55300121966SStefano Zampini       ierr = MatConjugate_SeqDense(A);CHKERRQ(ierr);
55400121966SStefano Zampini #endif
55500121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5568b83055fSJed Brown       PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info));
55700121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
55800121966SStefano Zampini #if defined(PETSC_USE_COMPLEX)
55900121966SStefano Zampini       ierr = MatConjugate_SeqDense(A);CHKERRQ(ierr);
56000121966SStefano Zampini #endif
561a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve");
562a49dc2a2SStefano Zampini #if defined(PETSC_USE_COMPLEX)
563a49dc2a2SStefano Zampini     } else if (A->hermitian) {
56400121966SStefano Zampini       ierr = MatConjugate_SeqDense(A);CHKERRQ(ierr);
56500121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
56600121966SStefano Zampini       PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
56700121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
56800121966SStefano Zampini       ierr = MatConjugate_SeqDense(A);CHKERRQ(ierr);
569ae7cfcebSSatish Balay #endif
570a49dc2a2SStefano Zampini     } else { /* symmetric case */
57100121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
572a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info));
57300121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
574a49dc2a2SStefano Zampini       if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve");
575da3a660dSBarry Smith     }
576a49dc2a2SStefano Zampini #endif
577a49dc2a2SStefano Zampini   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve");
578f1ceaac6SMatthew G. Knepley   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
5791ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
580dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
5813a40ed3dSBarry Smith   PetscFunctionReturn(0);
582da3a660dSBarry Smith }
5836ee01492SSatish Balay 
584e0877f53SBarry Smith static PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
585da3a660dSBarry Smith {
586dfbe8321SBarry Smith   PetscErrorCode    ierr;
587f1ceaac6SMatthew G. Knepley   const PetscScalar *x;
588f1ceaac6SMatthew G. Knepley   PetscScalar       *y,sone = 1.0;
589da3a660dSBarry Smith   Vec               tmp = 0;
59067e560aaSBarry Smith 
5913a40ed3dSBarry Smith   PetscFunctionBegin;
592f1ceaac6SMatthew G. Knepley   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
5931ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
594d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
595da3a660dSBarry Smith   if (yy == zz) {
59678b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
5973bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);CHKERRQ(ierr);
59878b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
599da3a660dSBarry Smith   }
6008208b9aeSStefano Zampini   ierr = MatSolve_SeqDense(A,xx,yy);CHKERRQ(ierr);
6016bf464f9SBarry Smith   if (tmp) {
6026bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
6036bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
6046bf464f9SBarry Smith   } else {
6056bf464f9SBarry Smith     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
6066bf464f9SBarry Smith   }
607f1ceaac6SMatthew G. Knepley   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
6081ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6098208b9aeSStefano Zampini   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
611da3a660dSBarry Smith }
61267e560aaSBarry Smith 
613e0877f53SBarry Smith static PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
614da3a660dSBarry Smith {
6156849ba73SBarry Smith   PetscErrorCode    ierr;
616f1ceaac6SMatthew G. Knepley   const PetscScalar *x;
617f1ceaac6SMatthew G. Knepley   PetscScalar       *y,sone = 1.0;
61889ae1891SBarry Smith   Vec               tmp = NULL;
61967e560aaSBarry Smith 
6203a40ed3dSBarry Smith   PetscFunctionBegin;
621d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
622f1ceaac6SMatthew G. Knepley   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6231ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
624da3a660dSBarry Smith   if (yy == zz) {
62578b31e54SBarry Smith     ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr);
6263bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);CHKERRQ(ierr);
62778b31e54SBarry Smith     ierr = VecCopy(yy,tmp);CHKERRQ(ierr);
628da3a660dSBarry Smith   }
6298208b9aeSStefano Zampini   ierr = MatSolveTranspose_SeqDense(A,xx,yy);CHKERRQ(ierr);
63090f02eecSBarry Smith   if (tmp) {
6312dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr);
6326bf464f9SBarry Smith     ierr = VecDestroy(&tmp);CHKERRQ(ierr);
6333a40ed3dSBarry Smith   } else {
6342dcb1b2aSMatthew Knepley     ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);
63590f02eecSBarry Smith   }
636f1ceaac6SMatthew G. Knepley   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
6371ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6388208b9aeSStefano Zampini   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
6393a40ed3dSBarry Smith   PetscFunctionReturn(0);
640da3a660dSBarry Smith }
641db4efbfdSBarry Smith 
642db4efbfdSBarry Smith /* ---------------------------------------------------------------*/
643db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots,
644db4efbfdSBarry Smith    rather than put it in the Mat->row slot.*/
645e0877f53SBarry Smith static PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo)
646db4efbfdSBarry Smith {
647db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF)
648db4efbfdSBarry Smith   PetscFunctionBegin;
649e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable.");
650db4efbfdSBarry Smith #else
651db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
652db4efbfdSBarry Smith   PetscErrorCode ierr;
653db4efbfdSBarry Smith   PetscBLASInt   n,m,info;
654db4efbfdSBarry Smith 
655db4efbfdSBarry Smith   PetscFunctionBegin;
656c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
657c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
658db4efbfdSBarry Smith   if (!mat->pivots) {
6598208b9aeSStefano Zampini     ierr = PetscMalloc1(A->rmap->n,&mat->pivots);CHKERRQ(ierr);
6603bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
661db4efbfdSBarry Smith   }
662db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
6638e57ea43SSatish Balay   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6648b83055fSJed Brown   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info));
6658e57ea43SSatish Balay   ierr = PetscFPTrapPop();CHKERRQ(ierr);
6668e57ea43SSatish Balay 
667e32f2f54SBarry Smith   if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization");
668e32f2f54SBarry Smith   if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization");
6698208b9aeSStefano Zampini 
670db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
6718208b9aeSStefano Zampini   A->ops->matsolve          = MatMatSolve_SeqDense;
672db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
673db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
674db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
675d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_LU;
676db4efbfdSBarry Smith 
677f6224b95SHong Zhang   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
678f6224b95SHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);CHKERRQ(ierr);
679f6224b95SHong Zhang 
680dc0b31edSSatish Balay   ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr);
681db4efbfdSBarry Smith #endif
682db4efbfdSBarry Smith   PetscFunctionReturn(0);
683db4efbfdSBarry Smith }
684db4efbfdSBarry Smith 
685a49dc2a2SStefano Zampini /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */
686e0877f53SBarry Smith static PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo)
687db4efbfdSBarry Smith {
688db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF)
689db4efbfdSBarry Smith   PetscFunctionBegin;
690e32f2f54SBarry Smith   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable.");
691db4efbfdSBarry Smith #else
692db4efbfdSBarry Smith   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
693db4efbfdSBarry Smith   PetscErrorCode ierr;
694c5df96a5SBarry Smith   PetscBLASInt   info,n;
695db4efbfdSBarry Smith 
696db4efbfdSBarry Smith   PetscFunctionBegin;
697c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
698db4efbfdSBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
699a49dc2a2SStefano Zampini   if (A->spd) {
70000121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7018b83055fSJed Brown     PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info));
70200121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
703a49dc2a2SStefano Zampini #if defined(PETSC_USE_COMPLEX)
704a49dc2a2SStefano Zampini   } else if (A->hermitian) {
705a49dc2a2SStefano Zampini     if (!mat->pivots) {
706a49dc2a2SStefano Zampini       ierr = PetscMalloc1(A->rmap->n,&mat->pivots);CHKERRQ(ierr);
707a49dc2a2SStefano Zampini       ierr = PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
708a49dc2a2SStefano Zampini     }
709a49dc2a2SStefano Zampini     if (!mat->fwork) {
710a49dc2a2SStefano Zampini       PetscScalar dummy;
711a49dc2a2SStefano Zampini 
712a49dc2a2SStefano Zampini       mat->lfwork = -1;
71300121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
714a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
71500121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
716a49dc2a2SStefano Zampini       mat->lfwork = (PetscInt)PetscRealPart(dummy);
717a49dc2a2SStefano Zampini       ierr = PetscMalloc1(mat->lfwork,&mat->fwork);CHKERRQ(ierr);
718a49dc2a2SStefano Zampini       ierr = PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));CHKERRQ(ierr);
719a49dc2a2SStefano Zampini     }
72000121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
721a49dc2a2SStefano Zampini     PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
72200121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
723a49dc2a2SStefano Zampini #endif
724a49dc2a2SStefano Zampini   } else { /* symmetric case */
725a49dc2a2SStefano Zampini     if (!mat->pivots) {
726a49dc2a2SStefano Zampini       ierr = PetscMalloc1(A->rmap->n,&mat->pivots);CHKERRQ(ierr);
727a49dc2a2SStefano Zampini       ierr = PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr);
728a49dc2a2SStefano Zampini     }
729a49dc2a2SStefano Zampini     if (!mat->fwork) {
730a49dc2a2SStefano Zampini       PetscScalar dummy;
731a49dc2a2SStefano Zampini 
732a49dc2a2SStefano Zampini       mat->lfwork = -1;
73300121966SStefano Zampini       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
734a49dc2a2SStefano Zampini       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info));
73500121966SStefano Zampini       ierr = PetscFPTrapPop();CHKERRQ(ierr);
736a49dc2a2SStefano Zampini       mat->lfwork = (PetscInt)PetscRealPart(dummy);
737a49dc2a2SStefano Zampini       ierr = PetscMalloc1(mat->lfwork,&mat->fwork);CHKERRQ(ierr);
738a49dc2a2SStefano Zampini       ierr = PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt));CHKERRQ(ierr);
739a49dc2a2SStefano Zampini     }
74000121966SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
741a49dc2a2SStefano Zampini     PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info));
74200121966SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
743a49dc2a2SStefano Zampini   }
744e32f2f54SBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1);
7458208b9aeSStefano Zampini 
746db4efbfdSBarry Smith   A->ops->solve             = MatSolve_SeqDense;
7478208b9aeSStefano Zampini   A->ops->matsolve          = MatMatSolve_SeqDense;
748db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqDense;
749db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqDense;
750db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
751d5f3da31SBarry Smith   A->factortype             = MAT_FACTOR_CHOLESKY;
7522205254eSKarl Rupp 
753f6224b95SHong Zhang   ierr = PetscFree(A->solvertype);CHKERRQ(ierr);
754f6224b95SHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&A->solvertype);CHKERRQ(ierr);
755f6224b95SHong Zhang 
756eb3f19e4SBarry Smith   ierr = PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr);
757db4efbfdSBarry Smith #endif
758db4efbfdSBarry Smith   PetscFunctionReturn(0);
759db4efbfdSBarry Smith }
760db4efbfdSBarry Smith 
761db4efbfdSBarry Smith 
7620481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy)
763db4efbfdSBarry Smith {
764db4efbfdSBarry Smith   PetscErrorCode ierr;
765db4efbfdSBarry Smith   MatFactorInfo  info;
766db4efbfdSBarry Smith 
767db4efbfdSBarry Smith   PetscFunctionBegin;
768db4efbfdSBarry Smith   info.fill = 1.0;
7692205254eSKarl Rupp 
770c3ef05f6SHong Zhang   ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr);
771719d5645SBarry Smith   ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr);
772db4efbfdSBarry Smith   PetscFunctionReturn(0);
773db4efbfdSBarry Smith }
774db4efbfdSBarry Smith 
775e0877f53SBarry Smith static PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info)
776db4efbfdSBarry Smith {
777db4efbfdSBarry Smith   PetscFunctionBegin;
778c3ef05f6SHong Zhang   fact->assembled                  = PETSC_TRUE;
7791bbcc794SSatish Balay   fact->preallocated               = PETSC_TRUE;
780719d5645SBarry Smith   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense;
781bd443b22SStefano Zampini   fact->ops->solve                 = MatSolve_SeqDense;
782bd443b22SStefano Zampini   fact->ops->matsolve              = MatMatSolve_SeqDense;
783bd443b22SStefano Zampini   fact->ops->solvetranspose        = MatSolveTranspose_SeqDense;
784bd443b22SStefano Zampini   fact->ops->solveadd              = MatSolveAdd_SeqDense;
785bd443b22SStefano Zampini   fact->ops->solvetransposeadd     = MatSolveTransposeAdd_SeqDense;
786db4efbfdSBarry Smith   PetscFunctionReturn(0);
787db4efbfdSBarry Smith }
788db4efbfdSBarry Smith 
789e0877f53SBarry Smith static PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
790db4efbfdSBarry Smith {
791db4efbfdSBarry Smith   PetscFunctionBegin;
792b66fe19dSMatthew G Knepley   fact->preallocated           = PETSC_TRUE;
793c3ef05f6SHong Zhang   fact->assembled              = PETSC_TRUE;
794719d5645SBarry Smith   fact->ops->lufactornumeric   = MatLUFactorNumeric_SeqDense;
795bd443b22SStefano Zampini   fact->ops->solve             = MatSolve_SeqDense;
796bd443b22SStefano Zampini   fact->ops->matsolve          = MatMatSolve_SeqDense;
797bd443b22SStefano Zampini   fact->ops->solvetranspose    = MatSolveTranspose_SeqDense;
798bd443b22SStefano Zampini   fact->ops->solveadd          = MatSolveAdd_SeqDense;
799bd443b22SStefano Zampini   fact->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense;
800db4efbfdSBarry Smith   PetscFunctionReturn(0);
801db4efbfdSBarry Smith }
802db4efbfdSBarry Smith 
803cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact)
804db4efbfdSBarry Smith {
805db4efbfdSBarry Smith   PetscErrorCode ierr;
806db4efbfdSBarry Smith 
807db4efbfdSBarry Smith   PetscFunctionBegin;
808ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),fact);CHKERRQ(ierr);
809db4efbfdSBarry Smith   ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
810db4efbfdSBarry Smith   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
811db4efbfdSBarry Smith   if (ftype == MAT_FACTOR_LU) {
812db4efbfdSBarry Smith     (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense;
813db4efbfdSBarry Smith   } else {
814db4efbfdSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense;
815db4efbfdSBarry Smith   }
816d5f3da31SBarry Smith   (*fact)->factortype = ftype;
81700c67f3bSHong Zhang 
81800c67f3bSHong Zhang   ierr = PetscFree((*fact)->solvertype);CHKERRQ(ierr);
81900c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype);CHKERRQ(ierr);
820db4efbfdSBarry Smith   PetscFunctionReturn(0);
821db4efbfdSBarry Smith }
822db4efbfdSBarry Smith 
823289bc588SBarry Smith /* ------------------------------------------------------------------*/
824e0877f53SBarry Smith static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx)
825289bc588SBarry Smith {
826c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
827d9ca1df4SBarry Smith   PetscScalar       *x,*v = mat->v,zero = 0.0,xt;
828d9ca1df4SBarry Smith   const PetscScalar *b;
829dfbe8321SBarry Smith   PetscErrorCode    ierr;
830d0f46423SBarry Smith   PetscInt          m = A->rmap->n,i;
831c5df96a5SBarry Smith   PetscBLASInt      o = 1,bm;
832289bc588SBarry Smith 
8333a40ed3dSBarry Smith   PetscFunctionBegin;
834422a814eSBarry Smith   if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */
835c5df96a5SBarry Smith   ierr = PetscBLASIntCast(m,&bm);CHKERRQ(ierr);
836289bc588SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
8373bffc371SBarry Smith     /* this is a hack fix, should have another version without the second BLASdotu */
8382dcb1b2aSMatthew Knepley     ierr = VecSet(xx,zero);CHKERRQ(ierr);
839289bc588SBarry Smith   }
8401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
841d9ca1df4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
842b965ef7fSBarry Smith   its  = its*lits;
843e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
844289bc588SBarry Smith   while (its--) {
845fccaa45eSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
846289bc588SBarry Smith       for (i=0; i<m; i++) {
8473bffc371SBarry Smith         PetscStackCallBLAS("BLASdotu",xt   = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
84855a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
849289bc588SBarry Smith       }
850289bc588SBarry Smith     }
851fccaa45eSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
852289bc588SBarry Smith       for (i=m-1; i>=0; i--) {
8533bffc371SBarry Smith         PetscStackCallBLAS("BLASdotu",xt   = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o));
85455a1b374SBarry Smith         x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift);
855289bc588SBarry Smith       }
856289bc588SBarry Smith     }
857289bc588SBarry Smith   }
858d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
8591ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861289bc588SBarry Smith }
862289bc588SBarry Smith 
863289bc588SBarry Smith /* -----------------------------------------------------------------*/
864cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy)
865289bc588SBarry Smith {
866c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
867d9ca1df4SBarry Smith   const PetscScalar *v   = mat->v,*x;
868d9ca1df4SBarry Smith   PetscScalar       *y;
869dfbe8321SBarry Smith   PetscErrorCode    ierr;
8700805154bSBarry Smith   PetscBLASInt      m, n,_One=1;
871ea709b57SSatish Balay   PetscScalar       _DOne=1.0,_DZero=0.0;
8723a40ed3dSBarry Smith 
8733a40ed3dSBarry Smith   PetscFunctionBegin;
874c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
875c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
876d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8771ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8785ac36cfcSBarry Smith   if (!A->rmap->n || !A->cmap->n) {
8795ac36cfcSBarry Smith     PetscBLASInt i;
8805ac36cfcSBarry Smith     for (i=0; i<n; i++) y[i] = 0.0;
8815ac36cfcSBarry Smith   } else {
8828b83055fSJed Brown     PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One));
8835ac36cfcSBarry Smith     ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr);
8845ac36cfcSBarry Smith   }
885d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8861ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8873a40ed3dSBarry Smith   PetscFunctionReturn(0);
888289bc588SBarry Smith }
889800995b7SMatthew Knepley 
890cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy)
891289bc588SBarry Smith {
892c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
893d9ca1df4SBarry Smith   PetscScalar       *y,_DOne=1.0,_DZero=0.0;
894dfbe8321SBarry Smith   PetscErrorCode    ierr;
8950805154bSBarry Smith   PetscBLASInt      m, n, _One=1;
896d9ca1df4SBarry Smith   const PetscScalar *v = mat->v,*x;
8973a40ed3dSBarry Smith 
8983a40ed3dSBarry Smith   PetscFunctionBegin;
899c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
900c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
901d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9021ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
9035ac36cfcSBarry Smith   if (!A->rmap->n || !A->cmap->n) {
9045ac36cfcSBarry Smith     PetscBLASInt i;
9055ac36cfcSBarry Smith     for (i=0; i<m; i++) y[i] = 0.0;
9065ac36cfcSBarry Smith   } else {
9078b83055fSJed Brown     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One));
9085ac36cfcSBarry Smith     ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr);
9095ac36cfcSBarry Smith   }
910d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9111ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9123a40ed3dSBarry Smith   PetscFunctionReturn(0);
913289bc588SBarry Smith }
9146ee01492SSatish Balay 
915cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
916289bc588SBarry Smith {
917c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
918d9ca1df4SBarry Smith   const PetscScalar *v = mat->v,*x;
919d9ca1df4SBarry Smith   PetscScalar       *y,_DOne=1.0;
920dfbe8321SBarry Smith   PetscErrorCode    ierr;
9210805154bSBarry Smith   PetscBLASInt      m, n, _One=1;
9223a40ed3dSBarry Smith 
9233a40ed3dSBarry Smith   PetscFunctionBegin;
924c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
925c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
926d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
927600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
928d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9291ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
9308b83055fSJed Brown   PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
931d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9321ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
933dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
9343a40ed3dSBarry Smith   PetscFunctionReturn(0);
935289bc588SBarry Smith }
9366ee01492SSatish Balay 
937e0877f53SBarry Smith static PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy)
938289bc588SBarry Smith {
939c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
940d9ca1df4SBarry Smith   const PetscScalar *v = mat->v,*x;
941d9ca1df4SBarry Smith   PetscScalar       *y;
942dfbe8321SBarry Smith   PetscErrorCode    ierr;
9430805154bSBarry Smith   PetscBLASInt      m, n, _One=1;
94487828ca2SBarry Smith   PetscScalar       _DOne=1.0;
9453a40ed3dSBarry Smith 
9463a40ed3dSBarry Smith   PetscFunctionBegin;
947c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr);
948c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr);
949d0f46423SBarry Smith   if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0);
950600d6d0bSBarry Smith   if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);}
951d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9521ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
9538b83055fSJed Brown   PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One));
954d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9551ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
956dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr);
9573a40ed3dSBarry Smith   PetscFunctionReturn(0);
958289bc588SBarry Smith }
959289bc588SBarry Smith 
960289bc588SBarry Smith /* -----------------------------------------------------------------*/
961e0877f53SBarry Smith static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
962289bc588SBarry Smith {
963c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
96487828ca2SBarry Smith   PetscScalar    *v;
9656849ba73SBarry Smith   PetscErrorCode ierr;
96613f74950SBarry Smith   PetscInt       i;
96767e560aaSBarry Smith 
9683a40ed3dSBarry Smith   PetscFunctionBegin;
969d0f46423SBarry Smith   *ncols = A->cmap->n;
970289bc588SBarry Smith   if (cols) {
971854ce69bSBarry Smith     ierr = PetscMalloc1(A->cmap->n+1,cols);CHKERRQ(ierr);
972d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) (*cols)[i] = i;
973289bc588SBarry Smith   }
974289bc588SBarry Smith   if (vals) {
975854ce69bSBarry Smith     ierr = PetscMalloc1(A->cmap->n+1,vals);CHKERRQ(ierr);
976289bc588SBarry Smith     v    = mat->v + row;
977d0f46423SBarry Smith     for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;}
978289bc588SBarry Smith   }
9793a40ed3dSBarry Smith   PetscFunctionReturn(0);
980289bc588SBarry Smith }
9816ee01492SSatish Balay 
982e0877f53SBarry Smith static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
983289bc588SBarry Smith {
984dfbe8321SBarry Smith   PetscErrorCode ierr;
9856e111a19SKarl Rupp 
986606d414cSSatish Balay   PetscFunctionBegin;
987606d414cSSatish Balay   if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);}
988606d414cSSatish Balay   if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); }
9893a40ed3dSBarry Smith   PetscFunctionReturn(0);
990289bc588SBarry Smith }
991289bc588SBarry Smith /* ----------------------------------------------------------------*/
992e0877f53SBarry Smith static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv)
993289bc588SBarry Smith {
994c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
99513f74950SBarry Smith   PetscInt     i,j,idx=0;
996d6dfbf8fSBarry Smith 
9973a40ed3dSBarry Smith   PetscFunctionBegin;
998289bc588SBarry Smith   if (!mat->roworiented) {
999dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
1000289bc588SBarry Smith       for (j=0; j<n; j++) {
1001cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
10022515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1003e32f2f54SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
100458804f6dSBarry Smith #endif
1005289bc588SBarry Smith         for (i=0; i<m; i++) {
1006cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
10072515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1008e32f2f54SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
100958804f6dSBarry Smith #endif
1010cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1011289bc588SBarry Smith         }
1012289bc588SBarry Smith       }
10133a40ed3dSBarry Smith     } else {
1014289bc588SBarry Smith       for (j=0; j<n; j++) {
1015cddbea37SSatish Balay         if (indexn[j] < 0) {idx += m; continue;}
10162515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1017e32f2f54SBarry Smith         if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
101858804f6dSBarry Smith #endif
1019289bc588SBarry Smith         for (i=0; i<m; i++) {
1020cddbea37SSatish Balay           if (indexm[i] < 0) {idx++; continue;}
10212515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1022e32f2f54SBarry Smith           if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
102358804f6dSBarry Smith #endif
1024cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1025289bc588SBarry Smith         }
1026289bc588SBarry Smith       }
1027289bc588SBarry Smith     }
10283a40ed3dSBarry Smith   } else {
1029dbb450caSBarry Smith     if (addv == INSERT_VALUES) {
1030e8d4e0b9SBarry Smith       for (i=0; i<m; i++) {
1031cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
10322515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1033e32f2f54SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
103458804f6dSBarry Smith #endif
1035e8d4e0b9SBarry Smith         for (j=0; j<n; j++) {
1036cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
10372515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1038e32f2f54SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
103958804f6dSBarry Smith #endif
1040cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++];
1041e8d4e0b9SBarry Smith         }
1042e8d4e0b9SBarry Smith       }
10433a40ed3dSBarry Smith     } else {
1044289bc588SBarry Smith       for (i=0; i<m; i++) {
1045cddbea37SSatish Balay         if (indexm[i] < 0) { idx += n; continue;}
10462515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1047e32f2f54SBarry Smith         if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1);
104858804f6dSBarry Smith #endif
1049289bc588SBarry Smith         for (j=0; j<n; j++) {
1050cddbea37SSatish Balay           if (indexn[j] < 0) { idx++; continue;}
10512515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
1052e32f2f54SBarry Smith           if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1);
105358804f6dSBarry Smith #endif
1054cddbea37SSatish Balay           mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++];
1055289bc588SBarry Smith         }
1056289bc588SBarry Smith       }
1057289bc588SBarry Smith     }
1058e8d4e0b9SBarry Smith   }
10593a40ed3dSBarry Smith   PetscFunctionReturn(0);
1060289bc588SBarry Smith }
1061e8d4e0b9SBarry Smith 
1062e0877f53SBarry Smith static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[])
1063ae80bb75SLois Curfman McInnes {
1064ae80bb75SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
106513f74950SBarry Smith   PetscInt     i,j;
1066ae80bb75SLois Curfman McInnes 
10673a40ed3dSBarry Smith   PetscFunctionBegin;
1068ae80bb75SLois Curfman McInnes   /* row-oriented output */
1069ae80bb75SLois Curfman McInnes   for (i=0; i<m; i++) {
107097e567efSBarry Smith     if (indexm[i] < 0) {v += n;continue;}
1071e32f2f54SBarry Smith     if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n);
1072ae80bb75SLois Curfman McInnes     for (j=0; j<n; j++) {
10736f31f424SBarry Smith       if (indexn[j] < 0) {v++; continue;}
1074e32f2f54SBarry Smith       if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n);
107597e567efSBarry Smith       *v++ = mat->v[indexn[j]*mat->lda + indexm[i]];
1076ae80bb75SLois Curfman McInnes     }
1077ae80bb75SLois Curfman McInnes   }
10783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1079ae80bb75SLois Curfman McInnes }
1080ae80bb75SLois Curfman McInnes 
1081289bc588SBarry Smith /* -----------------------------------------------------------------*/
1082289bc588SBarry Smith 
1083e0877f53SBarry Smith static PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer)
1084aabbc4fbSShri Abhyankar {
1085aabbc4fbSShri Abhyankar   Mat_SeqDense   *a;
1086aabbc4fbSShri Abhyankar   PetscErrorCode ierr;
1087aabbc4fbSShri Abhyankar   PetscInt       *scols,i,j,nz,header[4];
1088aabbc4fbSShri Abhyankar   int            fd;
1089aabbc4fbSShri Abhyankar   PetscMPIInt    size;
1090aabbc4fbSShri Abhyankar   PetscInt       *rowlengths = 0,M,N,*cols,grows,gcols;
1091aabbc4fbSShri Abhyankar   PetscScalar    *vals,*svals,*v,*w;
1092ce94432eSBarry Smith   MPI_Comm       comm;
10937f489da9SVaclav Hapla   PetscBool      isbinary;
1094aabbc4fbSShri Abhyankar 
1095aabbc4fbSShri Abhyankar   PetscFunctionBegin;
10967f489da9SVaclav Hapla   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
10977f489da9SVaclav Hapla   if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)newmat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newmat)->type_name);
10987f489da9SVaclav Hapla 
1099c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
1100c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
1101ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
1102aabbc4fbSShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1103aabbc4fbSShri Abhyankar   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
1104aabbc4fbSShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
11059860990eSLisandro Dalcin   ierr = PetscBinaryRead(fd,header,4,NULL,PETSC_INT);CHKERRQ(ierr);
1106aabbc4fbSShri Abhyankar   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
1107aabbc4fbSShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
1108aabbc4fbSShri Abhyankar 
1109aabbc4fbSShri Abhyankar   /* set global size if not set already*/
1110aabbc4fbSShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
1111aabbc4fbSShri Abhyankar     ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr);
1112aabbc4fbSShri Abhyankar   } else {
1113aabbc4fbSShri Abhyankar     /* if sizes and type are already set, check if the vector global sizes are correct */
1114aabbc4fbSShri Abhyankar     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
1115aabbc4fbSShri Abhyankar     if (M != grows ||  N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols);
1116aabbc4fbSShri Abhyankar   }
1117e6324fbbSBarry Smith   a = (Mat_SeqDense*)newmat->data;
1118e6324fbbSBarry Smith   if (!a->user_alloc) {
11190298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
1120e6324fbbSBarry Smith   }
1121aabbc4fbSShri Abhyankar 
1122aabbc4fbSShri Abhyankar   if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */
1123aabbc4fbSShri Abhyankar     a = (Mat_SeqDense*)newmat->data;
1124aabbc4fbSShri Abhyankar     v = a->v;
1125aabbc4fbSShri Abhyankar     /* Allocate some temp space to read in the values and then flip them
1126aabbc4fbSShri Abhyankar        from row major to column major */
1127854ce69bSBarry Smith     ierr = PetscMalloc1(M*N > 0 ? M*N : 1,&w);CHKERRQ(ierr);
1128aabbc4fbSShri Abhyankar     /* read in nonzero values */
11299860990eSLisandro Dalcin     ierr = PetscBinaryRead(fd,w,M*N,NULL,PETSC_SCALAR);CHKERRQ(ierr);
1130aabbc4fbSShri Abhyankar     /* now flip the values and store them in the matrix*/
1131aabbc4fbSShri Abhyankar     for (j=0; j<N; j++) {
1132aabbc4fbSShri Abhyankar       for (i=0; i<M; i++) {
1133aabbc4fbSShri Abhyankar         *v++ =w[i*N+j];
1134aabbc4fbSShri Abhyankar       }
1135aabbc4fbSShri Abhyankar     }
1136aabbc4fbSShri Abhyankar     ierr = PetscFree(w);CHKERRQ(ierr);
1137aabbc4fbSShri Abhyankar   } else {
1138aabbc4fbSShri Abhyankar     /* read row lengths */
1139854ce69bSBarry Smith     ierr = PetscMalloc1(M+1,&rowlengths);CHKERRQ(ierr);
11409860990eSLisandro Dalcin     ierr = PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);CHKERRQ(ierr);
1141aabbc4fbSShri Abhyankar 
1142aabbc4fbSShri Abhyankar     a = (Mat_SeqDense*)newmat->data;
1143aabbc4fbSShri Abhyankar     v = a->v;
1144aabbc4fbSShri Abhyankar 
1145aabbc4fbSShri Abhyankar     /* read column indices and nonzeros */
1146854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&scols);CHKERRQ(ierr);
1147aabbc4fbSShri Abhyankar     cols = scols;
11489860990eSLisandro Dalcin     ierr = PetscBinaryRead(fd,cols,nz,NULL,PETSC_INT);CHKERRQ(ierr);
1149854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&svals);CHKERRQ(ierr);
1150aabbc4fbSShri Abhyankar     vals = svals;
11519860990eSLisandro Dalcin     ierr = PetscBinaryRead(fd,vals,nz,NULL,PETSC_SCALAR);CHKERRQ(ierr);
1152aabbc4fbSShri Abhyankar 
1153aabbc4fbSShri Abhyankar     /* insert into matrix */
1154aabbc4fbSShri Abhyankar     for (i=0; i<M; i++) {
1155aabbc4fbSShri Abhyankar       for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j];
1156aabbc4fbSShri Abhyankar       svals += rowlengths[i]; scols += rowlengths[i];
1157aabbc4fbSShri Abhyankar     }
1158aabbc4fbSShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
1159aabbc4fbSShri Abhyankar     ierr = PetscFree(cols);CHKERRQ(ierr);
1160aabbc4fbSShri Abhyankar     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1161aabbc4fbSShri Abhyankar   }
1162aabbc4fbSShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1163aabbc4fbSShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1164aabbc4fbSShri Abhyankar   PetscFunctionReturn(0);
1165aabbc4fbSShri Abhyankar }
1166aabbc4fbSShri Abhyankar 
11676849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer)
1168289bc588SBarry Smith {
1169932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
1170dfbe8321SBarry Smith   PetscErrorCode    ierr;
117113f74950SBarry Smith   PetscInt          i,j;
11722dcb1b2aSMatthew Knepley   const char        *name;
117387828ca2SBarry Smith   PetscScalar       *v;
1174f3ef73ceSBarry Smith   PetscViewerFormat format;
11755f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX)
1176ace3abfcSBarry Smith   PetscBool         allreal = PETSC_TRUE;
11775f481a85SSatish Balay #endif
1178932b0c3eSLois Curfman McInnes 
11793a40ed3dSBarry Smith   PetscFunctionBegin;
1180b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1181456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
11823a40ed3dSBarry Smith     PetscFunctionReturn(0);  /* do nothing for now */
1183fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
1184d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1185d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
118644cd7ae7SLois Curfman McInnes       v    = a->v + i;
118777431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
1188d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1189aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1190329f5518SBarry Smith         if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) {
119157622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr);
1192329f5518SBarry Smith         } else if (PetscRealPart(*v)) {
119357622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));CHKERRQ(ierr);
11946831982aSBarry Smith         }
119580cd9d93SLois Curfman McInnes #else
11966831982aSBarry Smith         if (*v) {
119757622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);CHKERRQ(ierr);
11986831982aSBarry Smith         }
119980cd9d93SLois Curfman McInnes #endif
12001b807ce4Svictorle         v += a->lda;
120180cd9d93SLois Curfman McInnes       }
1202b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
120380cd9d93SLois Curfman McInnes     }
1204d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
12053a40ed3dSBarry Smith   } else {
1206d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
1207aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
120847989497SBarry Smith     /* determine if matrix has all real values */
120947989497SBarry Smith     v = a->v;
1210d0f46423SBarry Smith     for (i=0; i<A->rmap->n*A->cmap->n; i++) {
1211ffac6cdbSBarry Smith       if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;}
121247989497SBarry Smith     }
121347989497SBarry Smith #endif
1214fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
12153a7fca6bSBarry Smith       ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
1216d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1217d0f46423SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr);
1218fb9695e5SSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr);
1219ffac6cdbSBarry Smith     }
1220ffac6cdbSBarry Smith 
1221d0f46423SBarry Smith     for (i=0; i<A->rmap->n; i++) {
1222932b0c3eSLois Curfman McInnes       v = a->v + i;
1223d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1224aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
122547989497SBarry Smith         if (allreal) {
1226c61cd2faSJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));CHKERRQ(ierr);
122747989497SBarry Smith         } else {
1228c61cd2faSJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr);
122947989497SBarry Smith         }
1230289bc588SBarry Smith #else
1231c61cd2faSJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);CHKERRQ(ierr);
1232289bc588SBarry Smith #endif
12331b807ce4Svictorle         v += a->lda;
1234289bc588SBarry Smith       }
1235b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
1236289bc588SBarry Smith     }
1237fb9695e5SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1238b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr);
1239ffac6cdbSBarry Smith     }
1240d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
1241da3a660dSBarry Smith   }
1242b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
12433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1244289bc588SBarry Smith }
1245289bc588SBarry Smith 
12466849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer)
1247932b0c3eSLois Curfman McInnes {
1248932b0c3eSLois Curfman McInnes   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
12496849ba73SBarry Smith   PetscErrorCode    ierr;
125013f74950SBarry Smith   int               fd;
1251d0f46423SBarry Smith   PetscInt          ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n;
1252f4403165SShri Abhyankar   PetscScalar       *v,*anonz,*vals;
1253f4403165SShri Abhyankar   PetscViewerFormat format;
1254932b0c3eSLois Curfman McInnes 
12553a40ed3dSBarry Smith   PetscFunctionBegin;
1256b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
125790ace30eSBarry Smith 
1258f4403165SShri Abhyankar   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1259f4403165SShri Abhyankar   if (format == PETSC_VIEWER_NATIVE) {
1260f4403165SShri Abhyankar     /* store the matrix as a dense matrix */
1261785e854fSJed Brown     ierr = PetscMalloc1(4,&col_lens);CHKERRQ(ierr);
12622205254eSKarl Rupp 
1263f4403165SShri Abhyankar     col_lens[0] = MAT_FILE_CLASSID;
1264f4403165SShri Abhyankar     col_lens[1] = m;
1265f4403165SShri Abhyankar     col_lens[2] = n;
1266f4403165SShri Abhyankar     col_lens[3] = MATRIX_BINARY_FORMAT_DENSE;
12672205254eSKarl Rupp 
1268f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
1269f4403165SShri Abhyankar     ierr = PetscFree(col_lens);CHKERRQ(ierr);
1270f4403165SShri Abhyankar 
1271f4403165SShri Abhyankar     /* write out matrix, by rows */
1272854ce69bSBarry Smith     ierr = PetscMalloc1(m*n+1,&vals);CHKERRQ(ierr);
1273f4403165SShri Abhyankar     v    = a->v;
1274f4403165SShri Abhyankar     for (j=0; j<n; j++) {
1275f4403165SShri Abhyankar       for (i=0; i<m; i++) {
1276f4403165SShri Abhyankar         vals[j + i*n] = *v++;
1277f4403165SShri Abhyankar       }
1278f4403165SShri Abhyankar     }
1279f4403165SShri Abhyankar     ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
1280f4403165SShri Abhyankar     ierr = PetscFree(vals);CHKERRQ(ierr);
1281f4403165SShri Abhyankar   } else {
1282854ce69bSBarry Smith     ierr = PetscMalloc1(4+nz,&col_lens);CHKERRQ(ierr);
12832205254eSKarl Rupp 
12840700a824SBarry Smith     col_lens[0] = MAT_FILE_CLASSID;
1285932b0c3eSLois Curfman McInnes     col_lens[1] = m;
1286932b0c3eSLois Curfman McInnes     col_lens[2] = n;
1287932b0c3eSLois Curfman McInnes     col_lens[3] = nz;
1288932b0c3eSLois Curfman McInnes 
1289932b0c3eSLois Curfman McInnes     /* store lengths of each row and write (including header) to file */
1290932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) col_lens[4+i] = n;
12916f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr);
1292932b0c3eSLois Curfman McInnes 
1293932b0c3eSLois Curfman McInnes     /* Possibly should write in smaller increments, not whole matrix at once? */
1294932b0c3eSLois Curfman McInnes     /* store column indices (zero start index) */
1295932b0c3eSLois Curfman McInnes     ict = 0;
1296932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
1297932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) col_lens[ict++] = j;
1298932b0c3eSLois Curfman McInnes     }
12996f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1300606d414cSSatish Balay     ierr = PetscFree(col_lens);CHKERRQ(ierr);
1301932b0c3eSLois Curfman McInnes 
1302932b0c3eSLois Curfman McInnes     /* store nonzero values */
1303854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&anonz);CHKERRQ(ierr);
1304932b0c3eSLois Curfman McInnes     ict  = 0;
1305932b0c3eSLois Curfman McInnes     for (i=0; i<m; i++) {
1306932b0c3eSLois Curfman McInnes       v = a->v + i;
1307932b0c3eSLois Curfman McInnes       for (j=0; j<n; j++) {
13081b807ce4Svictorle         anonz[ict++] = *v; v += a->lda;
1309932b0c3eSLois Curfman McInnes       }
1310932b0c3eSLois Curfman McInnes     }
13116f69ff64SBarry Smith     ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr);
1312606d414cSSatish Balay     ierr = PetscFree(anonz);CHKERRQ(ierr);
1313f4403165SShri Abhyankar   }
13143a40ed3dSBarry Smith   PetscFunctionReturn(0);
1315932b0c3eSLois Curfman McInnes }
1316932b0c3eSLois Curfman McInnes 
13179804daf3SBarry Smith #include <petscdraw.h>
1318e0877f53SBarry Smith static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa)
1319f1af5d2fSBarry Smith {
1320f1af5d2fSBarry Smith   Mat               A  = (Mat) Aa;
1321f1af5d2fSBarry Smith   Mat_SeqDense      *a = (Mat_SeqDense*)A->data;
13226849ba73SBarry Smith   PetscErrorCode    ierr;
1323383922c3SLisandro Dalcin   PetscInt          m  = A->rmap->n,n = A->cmap->n,i,j;
1324383922c3SLisandro Dalcin   int               color = PETSC_DRAW_WHITE;
132587828ca2SBarry Smith   PetscScalar       *v = a->v;
1326b0a32e0cSBarry Smith   PetscViewer       viewer;
1327b05fc000SLisandro Dalcin   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r;
1328f3ef73ceSBarry Smith   PetscViewerFormat format;
1329f1af5d2fSBarry Smith 
1330f1af5d2fSBarry Smith   PetscFunctionBegin;
1331f1af5d2fSBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
1332b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1333b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1334f1af5d2fSBarry Smith 
1335f1af5d2fSBarry Smith   /* Loop over matrix elements drawing boxes */
1336383922c3SLisandro Dalcin 
1337fb9695e5SSatish Balay   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1338383922c3SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
1339f1af5d2fSBarry Smith     /* Blue for negative and Red for positive */
1340f1af5d2fSBarry Smith     for (j = 0; j < n; j++) {
1341383922c3SLisandro Dalcin       x_l = j; x_r = x_l + 1.0;
1342f1af5d2fSBarry Smith       for (i = 0; i < m; i++) {
1343f1af5d2fSBarry Smith         y_l = m - i - 1.0;
1344f1af5d2fSBarry Smith         y_r = y_l + 1.0;
1345329f5518SBarry Smith         if (PetscRealPart(v[j*m+i]) >  0.) {
1346b0a32e0cSBarry Smith           color = PETSC_DRAW_RED;
1347329f5518SBarry Smith         } else if (PetscRealPart(v[j*m+i]) <  0.) {
1348b0a32e0cSBarry Smith           color = PETSC_DRAW_BLUE;
1349f1af5d2fSBarry Smith         } else {
1350f1af5d2fSBarry Smith           continue;
1351f1af5d2fSBarry Smith         }
1352b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1353f1af5d2fSBarry Smith       }
1354f1af5d2fSBarry Smith     }
1355383922c3SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1356f1af5d2fSBarry Smith   } else {
1357f1af5d2fSBarry Smith     /* use contour shading to indicate magnitude of values */
1358f1af5d2fSBarry Smith     /* first determine max of all nonzero values */
1359b05fc000SLisandro Dalcin     PetscReal minv = 0.0, maxv = 0.0;
1360b05fc000SLisandro Dalcin     PetscDraw popup;
1361b05fc000SLisandro Dalcin 
1362f1af5d2fSBarry Smith     for (i=0; i < m*n; i++) {
1363f1af5d2fSBarry Smith       if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]);
1364f1af5d2fSBarry Smith     }
1365383922c3SLisandro Dalcin     if (minv >= maxv) maxv = minv + PETSC_SMALL;
1366b0a32e0cSBarry Smith     ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr);
136745f3bb6eSLisandro Dalcin     ierr = PetscDrawScalePopup(popup,minv,maxv);CHKERRQ(ierr);
1368383922c3SLisandro Dalcin 
1369383922c3SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
1370f1af5d2fSBarry Smith     for (j=0; j<n; j++) {
1371f1af5d2fSBarry Smith       x_l = j;
1372f1af5d2fSBarry Smith       x_r = x_l + 1.0;
1373f1af5d2fSBarry Smith       for (i=0; i<m; i++) {
1374f1af5d2fSBarry Smith         y_l = m - i - 1.0;
1375f1af5d2fSBarry Smith         y_r = y_l + 1.0;
1376b05fc000SLisandro Dalcin         color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv);
1377b0a32e0cSBarry Smith         ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr);
1378f1af5d2fSBarry Smith       }
1379f1af5d2fSBarry Smith     }
1380383922c3SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
1381f1af5d2fSBarry Smith   }
1382f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1383f1af5d2fSBarry Smith }
1384f1af5d2fSBarry Smith 
1385e0877f53SBarry Smith static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer)
1386f1af5d2fSBarry Smith {
1387b0a32e0cSBarry Smith   PetscDraw      draw;
1388ace3abfcSBarry Smith   PetscBool      isnull;
1389329f5518SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
1390dfbe8321SBarry Smith   PetscErrorCode ierr;
1391f1af5d2fSBarry Smith 
1392f1af5d2fSBarry Smith   PetscFunctionBegin;
1393b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1394b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
1395abc0a331SBarry Smith   if (isnull) PetscFunctionReturn(0);
1396f1af5d2fSBarry Smith 
1397d0f46423SBarry Smith   xr   = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
1398f1af5d2fSBarry Smith   xr  += w;          yr += h;        xl = -w;     yl = -h;
1399b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
1400832b7cebSLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
1401b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr);
14020298fd71SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr);
1403832b7cebSLisandro Dalcin   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
1404f1af5d2fSBarry Smith   PetscFunctionReturn(0);
1405f1af5d2fSBarry Smith }
1406f1af5d2fSBarry Smith 
1407dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer)
1408932b0c3eSLois Curfman McInnes {
1409dfbe8321SBarry Smith   PetscErrorCode ierr;
1410ace3abfcSBarry Smith   PetscBool      iascii,isbinary,isdraw;
1411932b0c3eSLois Curfman McInnes 
14123a40ed3dSBarry Smith   PetscFunctionBegin;
1413251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1414251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1415251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
14160f5bd95cSBarry Smith 
1417c45a1595SBarry Smith   if (iascii) {
1418c45a1595SBarry Smith     ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr);
14190f5bd95cSBarry Smith   } else if (isbinary) {
14203a40ed3dSBarry Smith     ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr);
1421f1af5d2fSBarry Smith   } else if (isdraw) {
1422f1af5d2fSBarry Smith     ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr);
1423932b0c3eSLois Curfman McInnes   }
14243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1425932b0c3eSLois Curfman McInnes }
1426289bc588SBarry Smith 
1427d3042a70SBarry Smith static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar array[])
1428d3042a70SBarry Smith {
1429d3042a70SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1430d3042a70SBarry Smith 
1431d3042a70SBarry Smith   PetscFunctionBegin;
1432d3042a70SBarry Smith   a->unplacedarray       = a->v;
1433d3042a70SBarry Smith   a->unplaced_user_alloc = a->user_alloc;
1434d3042a70SBarry Smith   a->v                   = (PetscScalar*) array;
1435d3042a70SBarry Smith   PetscFunctionReturn(0);
1436d3042a70SBarry Smith }
1437d3042a70SBarry Smith 
1438d3042a70SBarry Smith static PetscErrorCode MatDenseResetArray_SeqDense(Mat A)
1439d3042a70SBarry Smith {
1440d3042a70SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
1441d3042a70SBarry Smith 
1442d3042a70SBarry Smith   PetscFunctionBegin;
1443d3042a70SBarry Smith   a->v             = a->unplacedarray;
1444d3042a70SBarry Smith   a->user_alloc    = a->unplaced_user_alloc;
1445d3042a70SBarry Smith   a->unplacedarray = NULL;
1446d3042a70SBarry Smith   PetscFunctionReturn(0);
1447d3042a70SBarry Smith }
1448d3042a70SBarry Smith 
1449e0877f53SBarry Smith static PetscErrorCode MatDestroy_SeqDense(Mat mat)
1450289bc588SBarry Smith {
1451ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)mat->data;
1452dfbe8321SBarry Smith   PetscErrorCode ierr;
145390f02eecSBarry Smith 
14543a40ed3dSBarry Smith   PetscFunctionBegin;
1455aa482453SBarry Smith #if defined(PETSC_USE_LOG)
1456d0f46423SBarry Smith   PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n);
1457a5a9c739SBarry Smith #endif
145805b42c5fSBarry Smith   ierr = PetscFree(l->pivots);CHKERRQ(ierr);
1459a49dc2a2SStefano Zampini   ierr = PetscFree(l->fwork);CHKERRQ(ierr);
1460abc3b08eSStefano Zampini   ierr = MatDestroy(&l->ptapwork);CHKERRQ(ierr);
14616857c123SSatish Balay   if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);}
1462bf0cc555SLisandro Dalcin   ierr = PetscFree(mat->data);CHKERRQ(ierr);
1463dbd8c25aSHong Zhang 
1464dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
146549a6ff4bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL);CHKERRQ(ierr);
1466bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr);
1467d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL);CHKERRQ(ierr);
1468d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL);CHKERRQ(ierr);
1469bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr);
14708baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL);CHKERRQ(ierr);
14718baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
14728baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL);CHKERRQ(ierr);
14738baccfbdSHong Zhang #endif
1474bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);CHKERRQ(ierr);
1475bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
1476bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
1477bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
1478abc3b08eSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatPtAP_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
14793bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
14803bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
14813bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);CHKERRQ(ierr);
148286aefd0dSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL);CHKERRQ(ierr);
148386aefd0dSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL);CHKERRQ(ierr);
14843a40ed3dSBarry Smith   PetscFunctionReturn(0);
1485289bc588SBarry Smith }
1486289bc588SBarry Smith 
1487e0877f53SBarry Smith static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout)
1488289bc588SBarry Smith {
1489c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
14906849ba73SBarry Smith   PetscErrorCode ierr;
149113f74950SBarry Smith   PetscInt       k,j,m,n,M;
149287828ca2SBarry Smith   PetscScalar    *v,tmp;
149348b35521SBarry Smith 
14943a40ed3dSBarry Smith   PetscFunctionBegin;
1495d0f46423SBarry Smith   v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n;
14962847e3fdSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX && m == n) { /* in place transpose */
1497d3e5ee88SLois Curfman McInnes     for (j=0; j<m; j++) {
1498289bc588SBarry Smith       for (k=0; k<j; k++) {
14991b807ce4Svictorle         tmp        = v[j + k*M];
15001b807ce4Svictorle         v[j + k*M] = v[k + j*M];
15011b807ce4Svictorle         v[k + j*M] = tmp;
1502289bc588SBarry Smith       }
1503289bc588SBarry Smith     }
15043a40ed3dSBarry Smith   } else { /* out-of-place transpose */
1505d3e5ee88SLois Curfman McInnes     Mat          tmat;
1506ec8511deSBarry Smith     Mat_SeqDense *tmatd;
150787828ca2SBarry Smith     PetscScalar  *v2;
1508af36a384SStefano Zampini     PetscInt     M2;
1509ea709b57SSatish Balay 
15102847e3fdSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
1511ce94432eSBarry Smith       ierr = MatCreate(PetscObjectComm((PetscObject)A),&tmat);CHKERRQ(ierr);
1512d0f46423SBarry Smith       ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr);
15137adad957SLisandro Dalcin       ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
15140298fd71SBarry Smith       ierr = MatSeqDenseSetPreallocation(tmat,NULL);CHKERRQ(ierr);
1515fc4dec0aSBarry Smith     } else {
1516fc4dec0aSBarry Smith       tmat = *matout;
1517fc4dec0aSBarry Smith     }
1518ec8511deSBarry Smith     tmatd = (Mat_SeqDense*)tmat->data;
1519af36a384SStefano Zampini     v = mat->v; v2 = tmatd->v; M2 = tmatd->lda;
1520d3e5ee88SLois Curfman McInnes     for (j=0; j<n; j++) {
1521af36a384SStefano Zampini       for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M];
1522d3e5ee88SLois Curfman McInnes     }
15236d4a8577SBarry Smith     ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15246d4a8577SBarry Smith     ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15252847e3fdSStefano Zampini     if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) *matout = tmat;
15262847e3fdSStefano Zampini     else {
15272847e3fdSStefano Zampini       ierr = MatHeaderMerge(A,&tmat);CHKERRQ(ierr);
15282847e3fdSStefano Zampini     }
152948b35521SBarry Smith   }
15303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1531289bc588SBarry Smith }
1532289bc588SBarry Smith 
1533e0877f53SBarry Smith static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool  *flg)
1534289bc588SBarry Smith {
1535c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data;
1536c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data;
153713f74950SBarry Smith   PetscInt     i,j;
1538a2ea699eSBarry Smith   PetscScalar  *v1,*v2;
15399ea5d5aeSSatish Balay 
15403a40ed3dSBarry Smith   PetscFunctionBegin;
1541d0f46423SBarry Smith   if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1542d0f46423SBarry Smith   if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
1543d0f46423SBarry Smith   for (i=0; i<A1->rmap->n; i++) {
15441b807ce4Svictorle     v1 = mat1->v+i; v2 = mat2->v+i;
1545d0f46423SBarry Smith     for (j=0; j<A1->cmap->n; j++) {
15463a40ed3dSBarry Smith       if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);}
15471b807ce4Svictorle       v1 += mat1->lda; v2 += mat2->lda;
15481b807ce4Svictorle     }
1549289bc588SBarry Smith   }
155077c4ece6SBarry Smith   *flg = PETSC_TRUE;
15513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1552289bc588SBarry Smith }
1553289bc588SBarry Smith 
1554e0877f53SBarry Smith static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v)
1555289bc588SBarry Smith {
1556c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
1557dfbe8321SBarry Smith   PetscErrorCode ierr;
155813f74950SBarry Smith   PetscInt       i,n,len;
155987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
156044cd7ae7SLois Curfman McInnes 
15613a40ed3dSBarry Smith   PetscFunctionBegin;
15622dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
15637a97a34bSBarry Smith   ierr = VecGetSize(v,&n);CHKERRQ(ierr);
15641ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1565d0f46423SBarry Smith   len  = PetscMin(A->rmap->n,A->cmap->n);
1566e32f2f54SBarry Smith   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
156744cd7ae7SLois Curfman McInnes   for (i=0; i<len; i++) {
15681b807ce4Svictorle     x[i] = mat->v[i*mat->lda + i];
1569289bc588SBarry Smith   }
15701ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
15713a40ed3dSBarry Smith   PetscFunctionReturn(0);
1572289bc588SBarry Smith }
1573289bc588SBarry Smith 
1574e0877f53SBarry Smith static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr)
1575289bc588SBarry Smith {
1576c0bbcb79SLois Curfman McInnes   Mat_SeqDense      *mat = (Mat_SeqDense*)A->data;
1577f1ceaac6SMatthew G. Knepley   const PetscScalar *l,*r;
1578f1ceaac6SMatthew G. Knepley   PetscScalar       x,*v;
1579dfbe8321SBarry Smith   PetscErrorCode    ierr;
1580d0f46423SBarry Smith   PetscInt          i,j,m = A->rmap->n,n = A->cmap->n;
158155659b69SBarry Smith 
15823a40ed3dSBarry Smith   PetscFunctionBegin;
158328988994SBarry Smith   if (ll) {
15847a97a34bSBarry Smith     ierr = VecGetSize(ll,&m);CHKERRQ(ierr);
1585f1ceaac6SMatthew G. Knepley     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
1586e32f2f54SBarry Smith     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size");
1587da3a660dSBarry Smith     for (i=0; i<m; i++) {
1588da3a660dSBarry Smith       x = l[i];
1589da3a660dSBarry Smith       v = mat->v + i;
1590b43bac26SStefano Zampini       for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;}
1591da3a660dSBarry Smith     }
1592f1ceaac6SMatthew G. Knepley     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
1593eb3f19e4SBarry Smith     ierr = PetscLogFlops(1.0*n*m);CHKERRQ(ierr);
1594da3a660dSBarry Smith   }
159528988994SBarry Smith   if (rr) {
15967a97a34bSBarry Smith     ierr = VecGetSize(rr,&n);CHKERRQ(ierr);
1597f1ceaac6SMatthew G. Knepley     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
1598e32f2f54SBarry Smith     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size");
1599da3a660dSBarry Smith     for (i=0; i<n; i++) {
1600da3a660dSBarry Smith       x = r[i];
1601b43bac26SStefano Zampini       v = mat->v + i*mat->lda;
16022205254eSKarl Rupp       for (j=0; j<m; j++) (*v++) *= x;
1603da3a660dSBarry Smith     }
1604f1ceaac6SMatthew G. Knepley     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
1605eb3f19e4SBarry Smith     ierr = PetscLogFlops(1.0*n*m);CHKERRQ(ierr);
1606da3a660dSBarry Smith   }
16073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1608289bc588SBarry Smith }
1609289bc588SBarry Smith 
1610e0877f53SBarry Smith static PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm)
1611289bc588SBarry Smith {
1612c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
161387828ca2SBarry Smith   PetscScalar    *v   = mat->v;
1614329f5518SBarry Smith   PetscReal      sum  = 0.0;
1615d0f46423SBarry Smith   PetscInt       lda  =mat->lda,m=A->rmap->n,i,j;
1616efee365bSSatish Balay   PetscErrorCode ierr;
161755659b69SBarry Smith 
16183a40ed3dSBarry Smith   PetscFunctionBegin;
1619289bc588SBarry Smith   if (type == NORM_FROBENIUS) {
1620a5ce6ee0Svictorle     if (lda>m) {
1621d0f46423SBarry Smith       for (j=0; j<A->cmap->n; j++) {
1622a5ce6ee0Svictorle         v = mat->v+j*lda;
1623a5ce6ee0Svictorle         for (i=0; i<m; i++) {
1624a5ce6ee0Svictorle           sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1625a5ce6ee0Svictorle         }
1626a5ce6ee0Svictorle       }
1627a5ce6ee0Svictorle     } else {
1628570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FP16)
1629570b7f6dSBarry Smith       PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n;
1630570b7f6dSBarry Smith       *nrm = BLASnrm2_(&cnt,v,&one);
1631570b7f6dSBarry Smith     }
1632570b7f6dSBarry Smith #else
1633d0f46423SBarry Smith       for (i=0; i<A->cmap->n*A->rmap->n; i++) {
1634329f5518SBarry Smith         sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1635289bc588SBarry Smith       }
1636a5ce6ee0Svictorle     }
16378f1a2a5eSBarry Smith     *nrm = PetscSqrtReal(sum);
1638570b7f6dSBarry Smith #endif
1639dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
16403a40ed3dSBarry Smith   } else if (type == NORM_1) {
1641064f8208SBarry Smith     *nrm = 0.0;
1642d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
16431b807ce4Svictorle       v   = mat->v + j*mat->lda;
1644289bc588SBarry Smith       sum = 0.0;
1645d0f46423SBarry Smith       for (i=0; i<A->rmap->n; i++) {
164633a8263dSBarry Smith         sum += PetscAbsScalar(*v);  v++;
1647289bc588SBarry Smith       }
1648064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1649289bc588SBarry Smith     }
1650eb3f19e4SBarry Smith     ierr = PetscLogFlops(1.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
16513a40ed3dSBarry Smith   } else if (type == NORM_INFINITY) {
1652064f8208SBarry Smith     *nrm = 0.0;
1653d0f46423SBarry Smith     for (j=0; j<A->rmap->n; j++) {
1654289bc588SBarry Smith       v   = mat->v + j;
1655289bc588SBarry Smith       sum = 0.0;
1656d0f46423SBarry Smith       for (i=0; i<A->cmap->n; i++) {
16571b807ce4Svictorle         sum += PetscAbsScalar(*v); v += mat->lda;
1658289bc588SBarry Smith       }
1659064f8208SBarry Smith       if (sum > *nrm) *nrm = sum;
1660289bc588SBarry Smith     }
1661eb3f19e4SBarry Smith     ierr = PetscLogFlops(1.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr);
1662e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm");
16633a40ed3dSBarry Smith   PetscFunctionReturn(0);
1664289bc588SBarry Smith }
1665289bc588SBarry Smith 
1666e0877f53SBarry Smith static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg)
1667289bc588SBarry Smith {
1668c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *aij = (Mat_SeqDense*)A->data;
166963ba0a88SBarry Smith   PetscErrorCode ierr;
167067e560aaSBarry Smith 
16713a40ed3dSBarry Smith   PetscFunctionBegin;
1672b5a2b587SKris Buschelman   switch (op) {
1673b5a2b587SKris Buschelman   case MAT_ROW_ORIENTED:
16744e0d8c25SBarry Smith     aij->roworiented = flg;
1675b5a2b587SKris Buschelman     break;
1676512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1677b5a2b587SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
16783971808eSMatthew Knepley   case MAT_NEW_NONZERO_ALLOCATION_ERR:
16794e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
168013fa8e87SLisandro Dalcin   case MAT_KEEP_NONZERO_PATTERN:
1681b5a2b587SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1682b5a2b587SKris Buschelman   case MAT_USE_HASH_TABLE:
16830f8fb01aSBarry Smith   case MAT_IGNORE_ZERO_ENTRIES:
16845021d80fSJed Brown   case MAT_IGNORE_LOWER_TRIANGULAR:
1685*071fcb05SBarry Smith   case MAT_SORTED_FULL:
16865021d80fSJed Brown     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
16875021d80fSJed Brown     break;
16885021d80fSJed Brown   case MAT_SPD:
168977e54ba9SKris Buschelman   case MAT_SYMMETRIC:
169077e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
16919a4540c5SBarry Smith   case MAT_HERMITIAN:
16929a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
16935021d80fSJed Brown     /* These options are handled directly by MatSetOption() */
169477e54ba9SKris Buschelman     break;
1695b5a2b587SKris Buschelman   default:
1696e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]);
16973a40ed3dSBarry Smith   }
16983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1699289bc588SBarry Smith }
1700289bc588SBarry Smith 
1701e0877f53SBarry Smith static PetscErrorCode MatZeroEntries_SeqDense(Mat A)
17026f0a148fSBarry Smith {
1703ec8511deSBarry Smith   Mat_SeqDense   *l = (Mat_SeqDense*)A->data;
17046849ba73SBarry Smith   PetscErrorCode ierr;
1705d0f46423SBarry Smith   PetscInt       lda=l->lda,m=A->rmap->n,j;
17063a40ed3dSBarry Smith 
17073a40ed3dSBarry Smith   PetscFunctionBegin;
1708a5ce6ee0Svictorle   if (lda>m) {
1709d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
1710580bdb30SBarry Smith       ierr = PetscArrayzero(l->v+j*lda,m);CHKERRQ(ierr);
1711a5ce6ee0Svictorle     }
1712a5ce6ee0Svictorle   } else {
1713580bdb30SBarry Smith     ierr = PetscArrayzero(l->v,A->rmap->n*A->cmap->n);CHKERRQ(ierr);
1714a5ce6ee0Svictorle   }
17153a40ed3dSBarry Smith   PetscFunctionReturn(0);
17166f0a148fSBarry Smith }
17176f0a148fSBarry Smith 
1718e0877f53SBarry Smith static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
17196f0a148fSBarry Smith {
172097b48c8fSBarry Smith   PetscErrorCode    ierr;
1721ec8511deSBarry Smith   Mat_SeqDense      *l = (Mat_SeqDense*)A->data;
1722b9679d65SBarry Smith   PetscInt          m  = l->lda, n = A->cmap->n, i,j;
172397b48c8fSBarry Smith   PetscScalar       *slot,*bb;
172497b48c8fSBarry Smith   const PetscScalar *xx;
172555659b69SBarry Smith 
17263a40ed3dSBarry Smith   PetscFunctionBegin;
1727b9679d65SBarry Smith #if defined(PETSC_USE_DEBUG)
1728b9679d65SBarry Smith   for (i=0; i<N; i++) {
1729b9679d65SBarry Smith     if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed");
1730b9679d65SBarry Smith     if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n);
1731b9679d65SBarry Smith   }
1732b9679d65SBarry Smith #endif
1733b9679d65SBarry Smith 
173497b48c8fSBarry Smith   /* fix right hand side if needed */
173597b48c8fSBarry Smith   if (x && b) {
173697b48c8fSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
173797b48c8fSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
17382205254eSKarl Rupp     for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]];
173997b48c8fSBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
174097b48c8fSBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
174197b48c8fSBarry Smith   }
174297b48c8fSBarry Smith 
17436f0a148fSBarry Smith   for (i=0; i<N; i++) {
17446f0a148fSBarry Smith     slot = l->v + rows[i];
1745b9679d65SBarry Smith     for (j=0; j<n; j++) { *slot = 0.0; slot += m;}
17466f0a148fSBarry Smith   }
1747f4df32b1SMatthew Knepley   if (diag != 0.0) {
1748b9679d65SBarry Smith     if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices");
17496f0a148fSBarry Smith     for (i=0; i<N; i++) {
1750b9679d65SBarry Smith       slot  = l->v + (m+1)*rows[i];
1751f4df32b1SMatthew Knepley       *slot = diag;
17526f0a148fSBarry Smith     }
17536f0a148fSBarry Smith   }
17543a40ed3dSBarry Smith   PetscFunctionReturn(0);
17556f0a148fSBarry Smith }
1756557bce09SLois Curfman McInnes 
175749a6ff4bSBarry Smith static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda)
175849a6ff4bSBarry Smith {
175949a6ff4bSBarry Smith   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
176049a6ff4bSBarry Smith 
176149a6ff4bSBarry Smith   PetscFunctionBegin;
176249a6ff4bSBarry Smith   *lda = mat->lda;
176349a6ff4bSBarry Smith   PetscFunctionReturn(0);
176449a6ff4bSBarry Smith }
176549a6ff4bSBarry Smith 
1766e0877f53SBarry Smith static PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[])
176764e87e97SBarry Smith {
1768c0bbcb79SLois Curfman McInnes   Mat_SeqDense *mat = (Mat_SeqDense*)A->data;
17693a40ed3dSBarry Smith 
17703a40ed3dSBarry Smith   PetscFunctionBegin;
1771e32f2f54SBarry Smith   if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows");
177264e87e97SBarry Smith   *array = mat->v;
17733a40ed3dSBarry Smith   PetscFunctionReturn(0);
177464e87e97SBarry Smith }
17750754003eSLois Curfman McInnes 
1776e0877f53SBarry Smith static PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[])
1777ff14e315SSatish Balay {
17783a40ed3dSBarry Smith   PetscFunctionBegin;
17793a40ed3dSBarry Smith   PetscFunctionReturn(0);
1780ff14e315SSatish Balay }
17810754003eSLois Curfman McInnes 
1782dec5eb66SMatthew G Knepley /*@C
178349a6ff4bSBarry Smith    MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray()
178449a6ff4bSBarry Smith 
178549a6ff4bSBarry Smith    Logically Collective on Mat
178649a6ff4bSBarry Smith 
178749a6ff4bSBarry Smith    Input Parameter:
178849a6ff4bSBarry Smith .  mat - a MATSEQDENSE or MATMPIDENSE matrix
178949a6ff4bSBarry Smith 
179049a6ff4bSBarry Smith    Output Parameter:
179149a6ff4bSBarry Smith .   lda - the leading dimension
179249a6ff4bSBarry Smith 
179349a6ff4bSBarry Smith    Level: intermediate
179449a6ff4bSBarry Smith 
179549a6ff4bSBarry Smith .seealso: MatDenseGetArray(), MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead(), MatSeqDenseSetLDA()
179649a6ff4bSBarry Smith @*/
179749a6ff4bSBarry Smith PetscErrorCode  MatDenseGetLDA(Mat A,PetscInt *lda)
179849a6ff4bSBarry Smith {
179949a6ff4bSBarry Smith   PetscErrorCode ierr;
180049a6ff4bSBarry Smith 
180149a6ff4bSBarry Smith   PetscFunctionBegin;
180249a6ff4bSBarry Smith   ierr = PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda));CHKERRQ(ierr);
180349a6ff4bSBarry Smith   PetscFunctionReturn(0);
180449a6ff4bSBarry Smith }
180549a6ff4bSBarry Smith 
180649a6ff4bSBarry Smith /*@C
18078c778c55SBarry Smith    MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored
180873a71a0fSBarry Smith 
18098572280aSBarry Smith    Logically Collective on Mat
181073a71a0fSBarry Smith 
181173a71a0fSBarry Smith    Input Parameter:
1812579dbff0SBarry Smith .  mat - a MATSEQDENSE or MATMPIDENSE matrix
181373a71a0fSBarry Smith 
181473a71a0fSBarry Smith    Output Parameter:
181573a71a0fSBarry Smith .   array - pointer to the data
181673a71a0fSBarry Smith 
181773a71a0fSBarry Smith    Level: intermediate
181873a71a0fSBarry Smith 
18198572280aSBarry Smith .seealso: MatDenseRestoreArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
182073a71a0fSBarry Smith @*/
18218c778c55SBarry Smith PetscErrorCode  MatDenseGetArray(Mat A,PetscScalar **array)
182273a71a0fSBarry Smith {
182373a71a0fSBarry Smith   PetscErrorCode ierr;
182473a71a0fSBarry Smith 
182573a71a0fSBarry Smith   PetscFunctionBegin;
18268c778c55SBarry Smith   ierr = PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
182773a71a0fSBarry Smith   PetscFunctionReturn(0);
182873a71a0fSBarry Smith }
182973a71a0fSBarry Smith 
1830dec5eb66SMatthew G Knepley /*@C
1831579dbff0SBarry Smith    MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
183273a71a0fSBarry Smith 
18338572280aSBarry Smith    Logically Collective on Mat
18348572280aSBarry Smith 
18358572280aSBarry Smith    Input Parameters:
18368572280aSBarry Smith .  mat - a MATSEQDENSE or MATMPIDENSE matrix
18378572280aSBarry Smith .  array - pointer to the data
18388572280aSBarry Smith 
18398572280aSBarry Smith    Level: intermediate
18408572280aSBarry Smith 
18418572280aSBarry Smith .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArrayRead()
18428572280aSBarry Smith @*/
18438572280aSBarry Smith PetscErrorCode  MatDenseRestoreArray(Mat A,PetscScalar **array)
18448572280aSBarry Smith {
18458572280aSBarry Smith   PetscErrorCode ierr;
18468572280aSBarry Smith 
18478572280aSBarry Smith   PetscFunctionBegin;
18488572280aSBarry Smith   ierr = PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
18498572280aSBarry Smith   if (array) *array = NULL;
18508572280aSBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
18518572280aSBarry Smith   PetscFunctionReturn(0);
18528572280aSBarry Smith }
18538572280aSBarry Smith 
18548572280aSBarry Smith /*@C
18558572280aSBarry Smith    MatDenseGetArrayRead - gives access to the array where the data for a SeqDense matrix is stored
18568572280aSBarry Smith 
18578572280aSBarry Smith    Not Collective
18588572280aSBarry Smith 
18598572280aSBarry Smith    Input Parameter:
18608572280aSBarry Smith .  mat - a MATSEQDENSE or MATMPIDENSE matrix
18618572280aSBarry Smith 
18628572280aSBarry Smith    Output Parameter:
18638572280aSBarry Smith .   array - pointer to the data
18648572280aSBarry Smith 
18658572280aSBarry Smith    Level: intermediate
18668572280aSBarry Smith 
18678572280aSBarry Smith .seealso: MatDenseRestoreArray(), MatDenseGetArray(), MatDenseRestoreArrayRead()
18688572280aSBarry Smith @*/
18698572280aSBarry Smith PetscErrorCode  MatDenseGetArrayRead(Mat A,const PetscScalar **array)
18708572280aSBarry Smith {
18718572280aSBarry Smith   PetscErrorCode ierr;
18728572280aSBarry Smith 
18738572280aSBarry Smith   PetscFunctionBegin;
18748572280aSBarry Smith   ierr = PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array));CHKERRQ(ierr);
18758572280aSBarry Smith   PetscFunctionReturn(0);
18768572280aSBarry Smith }
18778572280aSBarry Smith 
18788572280aSBarry Smith /*@C
18798572280aSBarry Smith    MatDenseRestoreArrayRead - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray()
18808572280aSBarry Smith 
188173a71a0fSBarry Smith    Not Collective
188273a71a0fSBarry Smith 
188373a71a0fSBarry Smith    Input Parameters:
1884579dbff0SBarry Smith .  mat - a MATSEQDENSE or MATMPIDENSE matrix
188573a71a0fSBarry Smith .  array - pointer to the data
188673a71a0fSBarry Smith 
188773a71a0fSBarry Smith    Level: intermediate
188873a71a0fSBarry Smith 
18898572280aSBarry Smith .seealso: MatDenseGetArray(), MatDenseGetArrayRead(), MatDenseRestoreArray()
189073a71a0fSBarry Smith @*/
18918572280aSBarry Smith PetscErrorCode  MatDenseRestoreArrayRead(Mat A,const PetscScalar **array)
189273a71a0fSBarry Smith {
189373a71a0fSBarry Smith   PetscErrorCode ierr;
189473a71a0fSBarry Smith 
189573a71a0fSBarry Smith   PetscFunctionBegin;
18968572280aSBarry Smith   ierr = PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array));CHKERRQ(ierr);
18978572280aSBarry Smith   if (array) *array = NULL;
189873a71a0fSBarry Smith   PetscFunctionReturn(0);
189973a71a0fSBarry Smith }
190073a71a0fSBarry Smith 
19017dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
19020754003eSLois Curfman McInnes {
1903c0bbcb79SLois Curfman McInnes   Mat_SeqDense   *mat = (Mat_SeqDense*)A->data;
19046849ba73SBarry Smith   PetscErrorCode ierr;
19055d0c19d7SBarry Smith   PetscInt       i,j,nrows,ncols;
19065d0c19d7SBarry Smith   const PetscInt *irow,*icol;
190787828ca2SBarry Smith   PetscScalar    *av,*bv,*v = mat->v;
19080754003eSLois Curfman McInnes   Mat            newmat;
19090754003eSLois Curfman McInnes 
19103a40ed3dSBarry Smith   PetscFunctionBegin;
191178b31e54SBarry Smith   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
191278b31e54SBarry Smith   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
1913e03a110bSBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
1914e03a110bSBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
19150754003eSLois Curfman McInnes 
1916182d2002SSatish Balay   /* Check submatrixcall */
1917182d2002SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
191813f74950SBarry Smith     PetscInt n_cols,n_rows;
1919182d2002SSatish Balay     ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr);
192021a2c019SBarry Smith     if (n_rows != nrows || n_cols != ncols) {
1921f746d493SDmitry Karpeev       /* resize the result matrix to match number of requested rows/columns */
1922c61587bbSBarry Smith       ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
192321a2c019SBarry Smith     }
1924182d2002SSatish Balay     newmat = *B;
1925182d2002SSatish Balay   } else {
19260754003eSLois Curfman McInnes     /* Create and fill new matrix */
1927ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr);
1928f69a0ea3SMatthew Knepley     ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr);
19297adad957SLisandro Dalcin     ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr);
19300298fd71SBarry Smith     ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr);
1931182d2002SSatish Balay   }
1932182d2002SSatish Balay 
1933182d2002SSatish Balay   /* Now extract the data pointers and do the copy,column at a time */
1934182d2002SSatish Balay   bv = ((Mat_SeqDense*)newmat->data)->v;
1935182d2002SSatish Balay 
1936182d2002SSatish Balay   for (i=0; i<ncols; i++) {
19376de62eeeSBarry Smith     av = v + mat->lda*icol[i];
19382205254eSKarl Rupp     for (j=0; j<nrows; j++) *bv++ = av[irow[j]];
19390754003eSLois Curfman McInnes   }
1940182d2002SSatish Balay 
1941182d2002SSatish Balay   /* Assemble the matrices so that the correct flags are set */
19426d4a8577SBarry Smith   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19436d4a8577SBarry Smith   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19440754003eSLois Curfman McInnes 
19450754003eSLois Curfman McInnes   /* Free work space */
194678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
194778b31e54SBarry Smith   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
1948182d2002SSatish Balay   *B   = newmat;
19493a40ed3dSBarry Smith   PetscFunctionReturn(0);
19500754003eSLois Curfman McInnes }
19510754003eSLois Curfman McInnes 
19527dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1953905e6a2fSBarry Smith {
19546849ba73SBarry Smith   PetscErrorCode ierr;
195513f74950SBarry Smith   PetscInt       i;
1956905e6a2fSBarry Smith 
19573a40ed3dSBarry Smith   PetscFunctionBegin;
1958905e6a2fSBarry Smith   if (scall == MAT_INITIAL_MATRIX) {
1959df750dc8SHong Zhang     ierr = PetscCalloc1(n+1,B);CHKERRQ(ierr);
1960905e6a2fSBarry Smith   }
1961905e6a2fSBarry Smith 
1962905e6a2fSBarry Smith   for (i=0; i<n; i++) {
19637dae84e0SHong Zhang     ierr = MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
1964905e6a2fSBarry Smith   }
19653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1966905e6a2fSBarry Smith }
1967905e6a2fSBarry Smith 
1968e0877f53SBarry Smith static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode)
1969c0aa2d19SHong Zhang {
1970c0aa2d19SHong Zhang   PetscFunctionBegin;
1971c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1972c0aa2d19SHong Zhang }
1973c0aa2d19SHong Zhang 
1974e0877f53SBarry Smith static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode)
1975c0aa2d19SHong Zhang {
1976c0aa2d19SHong Zhang   PetscFunctionBegin;
1977c0aa2d19SHong Zhang   PetscFunctionReturn(0);
1978c0aa2d19SHong Zhang }
1979c0aa2d19SHong Zhang 
1980e0877f53SBarry Smith static PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str)
19814b0e389bSBarry Smith {
19824b0e389bSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data;
19836849ba73SBarry Smith   PetscErrorCode ierr;
1984d0f46423SBarry Smith   PetscInt       lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j;
19853a40ed3dSBarry Smith 
19863a40ed3dSBarry Smith   PetscFunctionBegin;
198733f4a19fSKris Buschelman   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
198833f4a19fSKris Buschelman   if (A->ops->copy != B->ops->copy) {
1989cb5b572fSBarry Smith     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
19903a40ed3dSBarry Smith     PetscFunctionReturn(0);
19913a40ed3dSBarry Smith   }
1992e32f2f54SBarry Smith   if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)");
1993a5ce6ee0Svictorle   if (lda1>m || lda2>m) {
19940dbb7854Svictorle     for (j=0; j<n; j++) {
1995580bdb30SBarry Smith       ierr = PetscArraycpy(b->v+j*lda2,a->v+j*lda1,m);CHKERRQ(ierr);
1996a5ce6ee0Svictorle     }
1997a5ce6ee0Svictorle   } else {
1998580bdb30SBarry Smith     ierr = PetscArraycpy(b->v,a->v,A->rmap->n*A->cmap->n);CHKERRQ(ierr);
1999a5ce6ee0Svictorle   }
2000cdc753b6SBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
2001273d9f13SBarry Smith   PetscFunctionReturn(0);
2002273d9f13SBarry Smith }
2003273d9f13SBarry Smith 
2004e0877f53SBarry Smith static PetscErrorCode MatSetUp_SeqDense(Mat A)
2005273d9f13SBarry Smith {
2006dfbe8321SBarry Smith   PetscErrorCode ierr;
2007273d9f13SBarry Smith 
2008273d9f13SBarry Smith   PetscFunctionBegin;
2009273d9f13SBarry Smith   ierr =  MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr);
20103a40ed3dSBarry Smith   PetscFunctionReturn(0);
20114b0e389bSBarry Smith }
20124b0e389bSBarry Smith 
2013ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A)
2014ba337c44SJed Brown {
2015ba337c44SJed Brown   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2016ba337c44SJed Brown   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2017ba337c44SJed Brown   PetscScalar  *aa = a->v;
2018ba337c44SJed Brown 
2019ba337c44SJed Brown   PetscFunctionBegin;
2020ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]);
2021ba337c44SJed Brown   PetscFunctionReturn(0);
2022ba337c44SJed Brown }
2023ba337c44SJed Brown 
2024ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A)
2025ba337c44SJed Brown {
2026ba337c44SJed Brown   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2027ba337c44SJed Brown   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2028ba337c44SJed Brown   PetscScalar  *aa = a->v;
2029ba337c44SJed Brown 
2030ba337c44SJed Brown   PetscFunctionBegin;
2031ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
2032ba337c44SJed Brown   PetscFunctionReturn(0);
2033ba337c44SJed Brown }
2034ba337c44SJed Brown 
2035ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A)
2036ba337c44SJed Brown {
2037ba337c44SJed Brown   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
2038ba337c44SJed Brown   PetscInt     i,nz = A->rmap->n*A->cmap->n;
2039ba337c44SJed Brown   PetscScalar  *aa = a->v;
2040ba337c44SJed Brown 
2041ba337c44SJed Brown   PetscFunctionBegin;
2042ba337c44SJed Brown   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
2043ba337c44SJed Brown   PetscFunctionReturn(0);
2044ba337c44SJed Brown }
2045284134d9SBarry Smith 
2046a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/
2047150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2048a9fe9ddaSSatish Balay {
2049a9fe9ddaSSatish Balay   PetscErrorCode ierr;
2050a9fe9ddaSSatish Balay 
2051a9fe9ddaSSatish Balay   PetscFunctionBegin;
2052a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
20533ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2054a9fe9ddaSSatish Balay     ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
20553ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2056a9fe9ddaSSatish Balay   }
20573ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
2058a9fe9ddaSSatish Balay   ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
20593ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
2060a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
2061a9fe9ddaSSatish Balay }
2062a9fe9ddaSSatish Balay 
2063a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2064a9fe9ddaSSatish Balay {
2065ee16a9a1SHong Zhang   PetscErrorCode ierr;
2066d0f46423SBarry Smith   PetscInt       m=A->rmap->n,n=B->cmap->n;
2067ee16a9a1SHong Zhang   Mat            Cmat;
2068a9fe9ddaSSatish Balay 
2069ee16a9a1SHong Zhang   PetscFunctionBegin;
2070e32f2f54SBarry Smith   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n);
2071ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
2072ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
2073ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
20740298fd71SBarry Smith   ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr);
2075d73949e8SHong Zhang 
2076ee16a9a1SHong Zhang   *C = Cmat;
2077ee16a9a1SHong Zhang   PetscFunctionReturn(0);
2078ee16a9a1SHong Zhang }
2079a9fe9ddaSSatish Balay 
2080a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2081a9fe9ddaSSatish Balay {
2082a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2083a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
2084a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
20850805154bSBarry Smith   PetscBLASInt   m,n,k;
2086a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
2087c5df96a5SBarry Smith   PetscErrorCode ierr;
2088fd4e9aacSBarry Smith   PetscBool      flg;
2089a9fe9ddaSSatish Balay 
2090a9fe9ddaSSatish Balay   PetscFunctionBegin;
2091fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr);
2092fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be dense");
2093fd4e9aacSBarry Smith 
2094fd4e9aacSBarry Smith   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
2095fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
2096fd4e9aacSBarry Smith   if (flg) {
2097fd4e9aacSBarry Smith     C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
2098fd4e9aacSBarry Smith     ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
2099fd4e9aacSBarry Smith     PetscFunctionReturn(0);
2100fd4e9aacSBarry Smith   }
2101fd4e9aacSBarry Smith 
2102fd4e9aacSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&flg);CHKERRQ(ierr);
2103fd4e9aacSBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be dense");
21048208b9aeSStefano Zampini   ierr = PetscBLASIntCast(C->rmap->n,&m);CHKERRQ(ierr);
21058208b9aeSStefano Zampini   ierr = PetscBLASIntCast(C->cmap->n,&n);CHKERRQ(ierr);
2106c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->cmap->n,&k);CHKERRQ(ierr);
210749d0e964SStefano Zampini   if (!m || !n || !k) PetscFunctionReturn(0);
21085ca1cc5dSStefano Zampini   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2109a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
2110a9fe9ddaSSatish Balay }
2111a9fe9ddaSSatish Balay 
211269f65d41SStefano Zampini PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
211369f65d41SStefano Zampini {
211469f65d41SStefano Zampini   PetscErrorCode ierr;
211569f65d41SStefano Zampini 
211669f65d41SStefano Zampini   PetscFunctionBegin;
211769f65d41SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
211869f65d41SStefano Zampini     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
211969f65d41SStefano Zampini     ierr = MatMatTransposeMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
212069f65d41SStefano Zampini     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
212169f65d41SStefano Zampini   }
212269f65d41SStefano Zampini   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
212369f65d41SStefano Zampini   ierr = MatMatTransposeMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
212469f65d41SStefano Zampini   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
212569f65d41SStefano Zampini   PetscFunctionReturn(0);
212669f65d41SStefano Zampini }
212769f65d41SStefano Zampini 
212869f65d41SStefano Zampini PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
212969f65d41SStefano Zampini {
213069f65d41SStefano Zampini   PetscErrorCode ierr;
213169f65d41SStefano Zampini   PetscInt       m=A->rmap->n,n=B->rmap->n;
213269f65d41SStefano Zampini   Mat            Cmat;
213369f65d41SStefano Zampini 
213469f65d41SStefano Zampini   PetscFunctionBegin;
213569f65d41SStefano Zampini   if (A->cmap->n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->cmap->n %d\n",A->cmap->n,B->cmap->n);
213669f65d41SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
213769f65d41SStefano Zampini   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
213869f65d41SStefano Zampini   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
213969f65d41SStefano Zampini   ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr);
214069f65d41SStefano Zampini 
214169f65d41SStefano Zampini   Cmat->assembled = PETSC_TRUE;
214269f65d41SStefano Zampini 
214369f65d41SStefano Zampini   *C = Cmat;
214469f65d41SStefano Zampini   PetscFunctionReturn(0);
214569f65d41SStefano Zampini }
214669f65d41SStefano Zampini 
214769f65d41SStefano Zampini PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
214869f65d41SStefano Zampini {
214969f65d41SStefano Zampini   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
215069f65d41SStefano Zampini   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
215169f65d41SStefano Zampini   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
215269f65d41SStefano Zampini   PetscBLASInt   m,n,k;
215369f65d41SStefano Zampini   PetscScalar    _DOne=1.0,_DZero=0.0;
215469f65d41SStefano Zampini   PetscErrorCode ierr;
215569f65d41SStefano Zampini 
215669f65d41SStefano Zampini   PetscFunctionBegin;
215749d0e964SStefano Zampini   ierr = PetscBLASIntCast(C->rmap->n,&m);CHKERRQ(ierr);
215849d0e964SStefano Zampini   ierr = PetscBLASIntCast(C->cmap->n,&n);CHKERRQ(ierr);
215969f65d41SStefano Zampini   ierr = PetscBLASIntCast(A->cmap->n,&k);CHKERRQ(ierr);
216049d0e964SStefano Zampini   if (!m || !n || !k) PetscFunctionReturn(0);
216169f65d41SStefano Zampini   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
216269f65d41SStefano Zampini   PetscFunctionReturn(0);
216369f65d41SStefano Zampini }
216469f65d41SStefano Zampini 
216575648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
2166a9fe9ddaSSatish Balay {
2167a9fe9ddaSSatish Balay   PetscErrorCode ierr;
2168a9fe9ddaSSatish Balay 
2169a9fe9ddaSSatish Balay   PetscFunctionBegin;
2170a9fe9ddaSSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21713ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
217275648e8dSHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr);
21733ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
2174a9fe9ddaSSatish Balay   }
21753ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
217675648e8dSHong Zhang   ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr);
21773ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
2178a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
2179a9fe9ddaSSatish Balay }
2180a9fe9ddaSSatish Balay 
218175648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
2182a9fe9ddaSSatish Balay {
2183ee16a9a1SHong Zhang   PetscErrorCode ierr;
2184d0f46423SBarry Smith   PetscInt       m=A->cmap->n,n=B->cmap->n;
2185ee16a9a1SHong Zhang   Mat            Cmat;
2186a9fe9ddaSSatish Balay 
2187ee16a9a1SHong Zhang   PetscFunctionBegin;
2188e32f2f54SBarry Smith   if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n);
2189ee16a9a1SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr);
2190ee16a9a1SHong Zhang   ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr);
2191ee16a9a1SHong Zhang   ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr);
21920298fd71SBarry Smith   ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr);
21932205254eSKarl Rupp 
2194ee16a9a1SHong Zhang   Cmat->assembled = PETSC_TRUE;
21952205254eSKarl Rupp 
2196ee16a9a1SHong Zhang   *C = Cmat;
2197ee16a9a1SHong Zhang   PetscFunctionReturn(0);
2198ee16a9a1SHong Zhang }
2199a9fe9ddaSSatish Balay 
220075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C)
2201a9fe9ddaSSatish Balay {
2202a9fe9ddaSSatish Balay   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2203a9fe9ddaSSatish Balay   Mat_SeqDense   *b = (Mat_SeqDense*)B->data;
2204a9fe9ddaSSatish Balay   Mat_SeqDense   *c = (Mat_SeqDense*)C->data;
22050805154bSBarry Smith   PetscBLASInt   m,n,k;
2206a9fe9ddaSSatish Balay   PetscScalar    _DOne=1.0,_DZero=0.0;
2207c5df96a5SBarry Smith   PetscErrorCode ierr;
2208a9fe9ddaSSatish Balay 
2209a9fe9ddaSSatish Balay   PetscFunctionBegin;
22108208b9aeSStefano Zampini   ierr = PetscBLASIntCast(C->rmap->n,&m);CHKERRQ(ierr);
22118208b9aeSStefano Zampini   ierr = PetscBLASIntCast(C->cmap->n,&n);CHKERRQ(ierr);
2212c5df96a5SBarry Smith   ierr = PetscBLASIntCast(A->rmap->n,&k);CHKERRQ(ierr);
221349d0e964SStefano Zampini   if (!m || !n || !k) PetscFunctionReturn(0);
22145ca1cc5dSStefano Zampini   PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda));
2215a9fe9ddaSSatish Balay   PetscFunctionReturn(0);
2216a9fe9ddaSSatish Balay }
2217985db425SBarry Smith 
2218e0877f53SBarry Smith static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[])
2219985db425SBarry Smith {
2220985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2221985db425SBarry Smith   PetscErrorCode ierr;
2222d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2223985db425SBarry Smith   PetscScalar    *x;
2224985db425SBarry Smith   MatScalar      *aa = a->v;
2225985db425SBarry Smith 
2226985db425SBarry Smith   PetscFunctionBegin;
2227e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2228985db425SBarry Smith 
2229985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
2230985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2231985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
2232e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2233985db425SBarry Smith   for (i=0; i<m; i++) {
2234985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
2235985db425SBarry Smith     for (j=1; j<n; j++) {
2236985db425SBarry Smith       if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
2237985db425SBarry Smith     }
2238985db425SBarry Smith   }
2239985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
2240985db425SBarry Smith   PetscFunctionReturn(0);
2241985db425SBarry Smith }
2242985db425SBarry Smith 
2243e0877f53SBarry Smith static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[])
2244985db425SBarry Smith {
2245985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2246985db425SBarry Smith   PetscErrorCode ierr;
2247d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2248985db425SBarry Smith   PetscScalar    *x;
2249985db425SBarry Smith   PetscReal      atmp;
2250985db425SBarry Smith   MatScalar      *aa = a->v;
2251985db425SBarry Smith 
2252985db425SBarry Smith   PetscFunctionBegin;
2253e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2254985db425SBarry Smith 
2255985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
2256985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2257985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
2258e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2259985db425SBarry Smith   for (i=0; i<m; i++) {
22609189402eSHong Zhang     x[i] = PetscAbsScalar(aa[i]);
2261985db425SBarry Smith     for (j=1; j<n; j++) {
2262985db425SBarry Smith       atmp = PetscAbsScalar(aa[i+m*j]);
2263985db425SBarry Smith       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;}
2264985db425SBarry Smith     }
2265985db425SBarry Smith   }
2266985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
2267985db425SBarry Smith   PetscFunctionReturn(0);
2268985db425SBarry Smith }
2269985db425SBarry Smith 
2270e0877f53SBarry Smith static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[])
2271985db425SBarry Smith {
2272985db425SBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
2273985db425SBarry Smith   PetscErrorCode ierr;
2274d0f46423SBarry Smith   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,p;
2275985db425SBarry Smith   PetscScalar    *x;
2276985db425SBarry Smith   MatScalar      *aa = a->v;
2277985db425SBarry Smith 
2278985db425SBarry Smith   PetscFunctionBegin;
2279e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2280985db425SBarry Smith 
2281985db425SBarry Smith   ierr = VecSet(v,0.0);CHKERRQ(ierr);
2282985db425SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2283985db425SBarry Smith   ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr);
2284e32f2f54SBarry Smith   if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2285985db425SBarry Smith   for (i=0; i<m; i++) {
2286985db425SBarry Smith     x[i] = aa[i]; if (idx) idx[i] = 0;
2287985db425SBarry Smith     for (j=1; j<n; j++) {
2288985db425SBarry Smith       if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;}
2289985db425SBarry Smith     }
2290985db425SBarry Smith   }
2291985db425SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
2292985db425SBarry Smith   PetscFunctionReturn(0);
2293985db425SBarry Smith }
2294985db425SBarry Smith 
2295e0877f53SBarry Smith static PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col)
22968d0534beSBarry Smith {
22978d0534beSBarry Smith   Mat_SeqDense   *a = (Mat_SeqDense*)A->data;
22988d0534beSBarry Smith   PetscErrorCode ierr;
22998d0534beSBarry Smith   PetscScalar    *x;
23008d0534beSBarry Smith 
23018d0534beSBarry Smith   PetscFunctionBegin;
2302e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
23038d0534beSBarry Smith 
23048d0534beSBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2305580bdb30SBarry Smith   ierr = PetscArraycpy(x,a->v+col*a->lda,A->rmap->n);CHKERRQ(ierr);
23068d0534beSBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
23078d0534beSBarry Smith   PetscFunctionReturn(0);
23088d0534beSBarry Smith }
23098d0534beSBarry Smith 
23100716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
23110716a85fSBarry Smith {
23120716a85fSBarry Smith   PetscErrorCode    ierr;
23130716a85fSBarry Smith   PetscInt          i,j,m,n;
23141683a169SBarry Smith   const PetscScalar *a;
23150716a85fSBarry Smith 
23160716a85fSBarry Smith   PetscFunctionBegin;
23170716a85fSBarry Smith   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
2318580bdb30SBarry Smith   ierr = PetscArrayzero(norms,n);CHKERRQ(ierr);
23191683a169SBarry Smith   ierr = MatDenseGetArrayRead(A,&a);CHKERRQ(ierr);
23200716a85fSBarry Smith   if (type == NORM_2) {
23210716a85fSBarry Smith     for (i=0; i<n; i++) {
23220716a85fSBarry Smith       for (j=0; j<m; j++) {
23230716a85fSBarry Smith         norms[i] += PetscAbsScalar(a[j]*a[j]);
23240716a85fSBarry Smith       }
23250716a85fSBarry Smith       a += m;
23260716a85fSBarry Smith     }
23270716a85fSBarry Smith   } else if (type == NORM_1) {
23280716a85fSBarry Smith     for (i=0; i<n; i++) {
23290716a85fSBarry Smith       for (j=0; j<m; j++) {
23300716a85fSBarry Smith         norms[i] += PetscAbsScalar(a[j]);
23310716a85fSBarry Smith       }
23320716a85fSBarry Smith       a += m;
23330716a85fSBarry Smith     }
23340716a85fSBarry Smith   } else if (type == NORM_INFINITY) {
23350716a85fSBarry Smith     for (i=0; i<n; i++) {
23360716a85fSBarry Smith       for (j=0; j<m; j++) {
23370716a85fSBarry Smith         norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
23380716a85fSBarry Smith       }
23390716a85fSBarry Smith       a += m;
23400716a85fSBarry Smith     }
2341ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType");
23421683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(A,&a);CHKERRQ(ierr);
23430716a85fSBarry Smith   if (type == NORM_2) {
23448f1a2a5eSBarry Smith     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
23450716a85fSBarry Smith   }
23460716a85fSBarry Smith   PetscFunctionReturn(0);
23470716a85fSBarry Smith }
23480716a85fSBarry Smith 
234973a71a0fSBarry Smith static PetscErrorCode  MatSetRandom_SeqDense(Mat x,PetscRandom rctx)
235073a71a0fSBarry Smith {
235173a71a0fSBarry Smith   PetscErrorCode ierr;
235273a71a0fSBarry Smith   PetscScalar    *a;
235373a71a0fSBarry Smith   PetscInt       m,n,i;
235473a71a0fSBarry Smith 
235573a71a0fSBarry Smith   PetscFunctionBegin;
235673a71a0fSBarry Smith   ierr = MatGetSize(x,&m,&n);CHKERRQ(ierr);
23578c778c55SBarry Smith   ierr = MatDenseGetArray(x,&a);CHKERRQ(ierr);
235873a71a0fSBarry Smith   for (i=0; i<m*n; i++) {
235973a71a0fSBarry Smith     ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr);
236073a71a0fSBarry Smith   }
23618c778c55SBarry Smith   ierr = MatDenseRestoreArray(x,&a);CHKERRQ(ierr);
236273a71a0fSBarry Smith   PetscFunctionReturn(0);
236373a71a0fSBarry Smith }
236473a71a0fSBarry Smith 
23653b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool  *missing,PetscInt *d)
23663b49f96aSBarry Smith {
23673b49f96aSBarry Smith   PetscFunctionBegin;
23683b49f96aSBarry Smith   *missing = PETSC_FALSE;
23693b49f96aSBarry Smith   PetscFunctionReturn(0);
23703b49f96aSBarry Smith }
237173a71a0fSBarry Smith 
2372af53bab2SHong Zhang static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals)
237386aefd0dSHong Zhang {
237486aefd0dSHong Zhang   Mat_SeqDense *a = (Mat_SeqDense*)A->data;
237586aefd0dSHong Zhang 
237686aefd0dSHong Zhang   PetscFunctionBegin;
237786aefd0dSHong Zhang   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
237886aefd0dSHong Zhang   *vals = a->v+col*a->lda;
237986aefd0dSHong Zhang   PetscFunctionReturn(0);
238086aefd0dSHong Zhang }
238186aefd0dSHong Zhang 
2382af53bab2SHong Zhang static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals)
238386aefd0dSHong Zhang {
238486aefd0dSHong Zhang   PetscFunctionBegin;
238586aefd0dSHong Zhang   *vals = 0; /* user cannot accidently use the array later */
238686aefd0dSHong Zhang   PetscFunctionReturn(0);
238786aefd0dSHong Zhang }
2388abc3b08eSStefano Zampini 
2389289bc588SBarry Smith /* -------------------------------------------------------------------*/
2390a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_SeqDense,
2391905e6a2fSBarry Smith                                         MatGetRow_SeqDense,
2392905e6a2fSBarry Smith                                         MatRestoreRow_SeqDense,
2393905e6a2fSBarry Smith                                         MatMult_SeqDense,
239497304618SKris Buschelman                                 /*  4*/ MatMultAdd_SeqDense,
23957c922b88SBarry Smith                                         MatMultTranspose_SeqDense,
23967c922b88SBarry Smith                                         MatMultTransposeAdd_SeqDense,
2397db4efbfdSBarry Smith                                         0,
2398db4efbfdSBarry Smith                                         0,
2399db4efbfdSBarry Smith                                         0,
2400db4efbfdSBarry Smith                                 /* 10*/ 0,
2401905e6a2fSBarry Smith                                         MatLUFactor_SeqDense,
2402905e6a2fSBarry Smith                                         MatCholeskyFactor_SeqDense,
240341f059aeSBarry Smith                                         MatSOR_SeqDense,
2404ec8511deSBarry Smith                                         MatTranspose_SeqDense,
240597304618SKris Buschelman                                 /* 15*/ MatGetInfo_SeqDense,
2406905e6a2fSBarry Smith                                         MatEqual_SeqDense,
2407905e6a2fSBarry Smith                                         MatGetDiagonal_SeqDense,
2408905e6a2fSBarry Smith                                         MatDiagonalScale_SeqDense,
2409905e6a2fSBarry Smith                                         MatNorm_SeqDense,
2410c0aa2d19SHong Zhang                                 /* 20*/ MatAssemblyBegin_SeqDense,
2411c0aa2d19SHong Zhang                                         MatAssemblyEnd_SeqDense,
2412905e6a2fSBarry Smith                                         MatSetOption_SeqDense,
2413905e6a2fSBarry Smith                                         MatZeroEntries_SeqDense,
2414d519adbfSMatthew Knepley                                 /* 24*/ MatZeroRows_SeqDense,
2415db4efbfdSBarry Smith                                         0,
2416db4efbfdSBarry Smith                                         0,
2417db4efbfdSBarry Smith                                         0,
2418db4efbfdSBarry Smith                                         0,
24194994cf47SJed Brown                                 /* 29*/ MatSetUp_SeqDense,
2420273d9f13SBarry Smith                                         0,
2421905e6a2fSBarry Smith                                         0,
242273a71a0fSBarry Smith                                         0,
242373a71a0fSBarry Smith                                         0,
2424d519adbfSMatthew Knepley                                 /* 34*/ MatDuplicate_SeqDense,
2425a5ae1ecdSBarry Smith                                         0,
2426a5ae1ecdSBarry Smith                                         0,
2427a5ae1ecdSBarry Smith                                         0,
2428a5ae1ecdSBarry Smith                                         0,
2429d519adbfSMatthew Knepley                                 /* 39*/ MatAXPY_SeqDense,
24307dae84e0SHong Zhang                                         MatCreateSubMatrices_SeqDense,
2431a5ae1ecdSBarry Smith                                         0,
24324b0e389bSBarry Smith                                         MatGetValues_SeqDense,
2433a5ae1ecdSBarry Smith                                         MatCopy_SeqDense,
2434d519adbfSMatthew Knepley                                 /* 44*/ MatGetRowMax_SeqDense,
2435a5ae1ecdSBarry Smith                                         MatScale_SeqDense,
24367d68702bSBarry Smith                                         MatShift_Basic,
2437a5ae1ecdSBarry Smith                                         0,
24383f49a652SStefano Zampini                                         MatZeroRowsColumns_SeqDense,
243973a71a0fSBarry Smith                                 /* 49*/ MatSetRandom_SeqDense,
2440a5ae1ecdSBarry Smith                                         0,
2441a5ae1ecdSBarry Smith                                         0,
2442a5ae1ecdSBarry Smith                                         0,
2443a5ae1ecdSBarry Smith                                         0,
2444d519adbfSMatthew Knepley                                 /* 54*/ 0,
2445a5ae1ecdSBarry Smith                                         0,
2446a5ae1ecdSBarry Smith                                         0,
2447a5ae1ecdSBarry Smith                                         0,
2448a5ae1ecdSBarry Smith                                         0,
2449d519adbfSMatthew Knepley                                 /* 59*/ 0,
2450e03a110bSBarry Smith                                         MatDestroy_SeqDense,
2451e03a110bSBarry Smith                                         MatView_SeqDense,
2452357abbc8SBarry Smith                                         0,
245397304618SKris Buschelman                                         0,
2454d519adbfSMatthew Knepley                                 /* 64*/ 0,
245597304618SKris Buschelman                                         0,
245697304618SKris Buschelman                                         0,
245797304618SKris Buschelman                                         0,
245897304618SKris Buschelman                                         0,
2459d519adbfSMatthew Knepley                                 /* 69*/ MatGetRowMaxAbs_SeqDense,
246097304618SKris Buschelman                                         0,
246197304618SKris Buschelman                                         0,
246297304618SKris Buschelman                                         0,
246397304618SKris Buschelman                                         0,
2464d519adbfSMatthew Knepley                                 /* 74*/ 0,
246597304618SKris Buschelman                                         0,
246697304618SKris Buschelman                                         0,
246797304618SKris Buschelman                                         0,
246897304618SKris Buschelman                                         0,
2469d519adbfSMatthew Knepley                                 /* 79*/ 0,
247097304618SKris Buschelman                                         0,
247197304618SKris Buschelman                                         0,
247297304618SKris Buschelman                                         0,
24735bba2384SShri Abhyankar                                 /* 83*/ MatLoad_SeqDense,
2474865e5f61SKris Buschelman                                         0,
24751cbb95d3SBarry Smith                                         MatIsHermitian_SeqDense,
2476865e5f61SKris Buschelman                                         0,
2477865e5f61SKris Buschelman                                         0,
2478865e5f61SKris Buschelman                                         0,
2479d519adbfSMatthew Knepley                                 /* 89*/ MatMatMult_SeqDense_SeqDense,
2480a9fe9ddaSSatish Balay                                         MatMatMultSymbolic_SeqDense_SeqDense,
2481a9fe9ddaSSatish Balay                                         MatMatMultNumeric_SeqDense_SeqDense,
2482abc3b08eSStefano Zampini                                         MatPtAP_SeqDense_SeqDense,
2483abc3b08eSStefano Zampini                                         MatPtAPSymbolic_SeqDense_SeqDense,
2484abc3b08eSStefano Zampini                                 /* 94*/ MatPtAPNumeric_SeqDense_SeqDense,
248569f65d41SStefano Zampini                                         MatMatTransposeMult_SeqDense_SeqDense,
248669f65d41SStefano Zampini                                         MatMatTransposeMultSymbolic_SeqDense_SeqDense,
248769f65d41SStefano Zampini                                         MatMatTransposeMultNumeric_SeqDense_SeqDense,
2488284134d9SBarry Smith                                         0,
2489d519adbfSMatthew Knepley                                 /* 99*/ 0,
2490284134d9SBarry Smith                                         0,
2491284134d9SBarry Smith                                         0,
2492ba337c44SJed Brown                                         MatConjugate_SeqDense,
2493f73d5cc4SBarry Smith                                         0,
2494ba337c44SJed Brown                                 /*104*/ 0,
2495ba337c44SJed Brown                                         MatRealPart_SeqDense,
2496ba337c44SJed Brown                                         MatImaginaryPart_SeqDense,
2497985db425SBarry Smith                                         0,
2498985db425SBarry Smith                                         0,
24998208b9aeSStefano Zampini                                 /*109*/ 0,
2500985db425SBarry Smith                                         0,
25018d0534beSBarry Smith                                         MatGetRowMin_SeqDense,
2502aabbc4fbSShri Abhyankar                                         MatGetColumnVector_SeqDense,
25033b49f96aSBarry Smith                                         MatMissingDiagonal_SeqDense,
2504aabbc4fbSShri Abhyankar                                 /*114*/ 0,
2505aabbc4fbSShri Abhyankar                                         0,
2506aabbc4fbSShri Abhyankar                                         0,
2507aabbc4fbSShri Abhyankar                                         0,
2508aabbc4fbSShri Abhyankar                                         0,
2509aabbc4fbSShri Abhyankar                                 /*119*/ 0,
2510aabbc4fbSShri Abhyankar                                         0,
2511aabbc4fbSShri Abhyankar                                         0,
25120716a85fSBarry Smith                                         0,
25130716a85fSBarry Smith                                         0,
25140716a85fSBarry Smith                                 /*124*/ 0,
25155df89d91SHong Zhang                                         MatGetColumnNorms_SeqDense,
25165df89d91SHong Zhang                                         0,
25175df89d91SHong Zhang                                         0,
25185df89d91SHong Zhang                                         0,
25195df89d91SHong Zhang                                 /*129*/ 0,
252075648e8dSHong Zhang                                         MatTransposeMatMult_SeqDense_SeqDense,
252175648e8dSHong Zhang                                         MatTransposeMatMultSymbolic_SeqDense_SeqDense,
252275648e8dSHong Zhang                                         MatTransposeMatMultNumeric_SeqDense_SeqDense,
25233964eb88SJed Brown                                         0,
25243964eb88SJed Brown                                 /*134*/ 0,
25253964eb88SJed Brown                                         0,
25263964eb88SJed Brown                                         0,
25273964eb88SJed Brown                                         0,
25283964eb88SJed Brown                                         0,
25293964eb88SJed Brown                                 /*139*/ 0,
2530f9426fe0SMark Adams                                         0,
2531d528f656SJakub Kruzik                                         0,
2532d528f656SJakub Kruzik                                         0,
2533d528f656SJakub Kruzik                                         0,
2534d528f656SJakub Kruzik                                 /*144*/ MatCreateMPIMatConcatenateSeqMat_SeqDense
2535985db425SBarry Smith };
253690ace30eSBarry Smith 
25374b828684SBarry Smith /*@C
2538fafbff53SBarry Smith    MatCreateSeqDense - Creates a sequential dense matrix that
2539d65003e9SLois Curfman McInnes    is stored in column major order (the usual Fortran 77 manner). Many
2540d65003e9SLois Curfman McInnes    of the matrix operations use the BLAS and LAPACK routines.
2541289bc588SBarry Smith 
2542d083f849SBarry Smith    Collective
2543db81eaa0SLois Curfman McInnes 
254420563c6bSBarry Smith    Input Parameters:
2545db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
25460c775827SLois Curfman McInnes .  m - number of rows
254718f449edSLois Curfman McInnes .  n - number of columns
25480298fd71SBarry Smith -  data - optional location of matrix data in column major order.  Set data=NULL for PETSc
2549dfc5480cSLois Curfman McInnes    to control all matrix memory allocation.
255020563c6bSBarry Smith 
255120563c6bSBarry Smith    Output Parameter:
255244cd7ae7SLois Curfman McInnes .  A - the matrix
255320563c6bSBarry Smith 
2554b259b22eSLois Curfman McInnes    Notes:
255518f449edSLois Curfman McInnes    The data input variable is intended primarily for Fortran programmers
255618f449edSLois Curfman McInnes    who wish to allocate their own matrix memory space.  Most users should
25570298fd71SBarry Smith    set data=NULL.
255818f449edSLois Curfman McInnes 
2559027ccd11SLois Curfman McInnes    Level: intermediate
2560027ccd11SLois Curfman McInnes 
256169b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues()
256220563c6bSBarry Smith @*/
25637087cfbeSBarry Smith PetscErrorCode  MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A)
2564289bc588SBarry Smith {
2565dfbe8321SBarry Smith   PetscErrorCode ierr;
25663b2fbd54SBarry Smith 
25673a40ed3dSBarry Smith   PetscFunctionBegin;
2568f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2569f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2570273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr);
2571273d9f13SBarry Smith   ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr);
2572273d9f13SBarry Smith   PetscFunctionReturn(0);
2573273d9f13SBarry Smith }
2574273d9f13SBarry Smith 
2575273d9f13SBarry Smith /*@C
2576273d9f13SBarry Smith    MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements
2577273d9f13SBarry Smith 
2578d083f849SBarry Smith    Collective
2579273d9f13SBarry Smith 
2580273d9f13SBarry Smith    Input Parameters:
25811c4f3114SJed Brown +  B - the matrix
25820298fd71SBarry Smith -  data - the array (or NULL)
2583273d9f13SBarry Smith 
2584273d9f13SBarry Smith    Notes:
2585273d9f13SBarry Smith    The data input variable is intended primarily for Fortran programmers
2586273d9f13SBarry Smith    who wish to allocate their own matrix memory space.  Most users should
2587284134d9SBarry Smith    need not call this routine.
2588273d9f13SBarry Smith 
2589273d9f13SBarry Smith    Level: intermediate
2590273d9f13SBarry Smith 
259169b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA()
2592867c911aSBarry Smith 
2593273d9f13SBarry Smith @*/
25947087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation(Mat B,PetscScalar data[])
2595273d9f13SBarry Smith {
25964ac538c5SBarry Smith   PetscErrorCode ierr;
2597a23d5eceSKris Buschelman 
2598a23d5eceSKris Buschelman   PetscFunctionBegin;
25994ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr);
2600a23d5eceSKris Buschelman   PetscFunctionReturn(0);
2601a23d5eceSKris Buschelman }
2602a23d5eceSKris Buschelman 
26037087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data)
2604a23d5eceSKris Buschelman {
2605273d9f13SBarry Smith   Mat_SeqDense   *b;
2606dfbe8321SBarry Smith   PetscErrorCode ierr;
2607273d9f13SBarry Smith 
2608273d9f13SBarry Smith   PetscFunctionBegin;
2609273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
2610a868139aSShri Abhyankar 
261134ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
261234ef9618SShri Abhyankar   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
261334ef9618SShri Abhyankar 
2614273d9f13SBarry Smith   b       = (Mat_SeqDense*)B->data;
261586d161a7SShri Abhyankar   b->Mmax = B->rmap->n;
261686d161a7SShri Abhyankar   b->Nmax = B->cmap->n;
261786d161a7SShri Abhyankar   if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n;
261886d161a7SShri Abhyankar 
2619220afb94SBarry Smith   ierr = PetscIntMultError(b->lda,b->Nmax,NULL);CHKERRQ(ierr);
26209e8f95c4SLisandro Dalcin   if (!data) { /* petsc-allocated storage */
26219e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2622e92229d0SSatish Balay     ierr = PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);CHKERRQ(ierr);
26233bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr);
26242205254eSKarl Rupp 
26259e8f95c4SLisandro Dalcin     b->user_alloc = PETSC_FALSE;
2626273d9f13SBarry Smith   } else { /* user-allocated storage */
26279e8f95c4SLisandro Dalcin     if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); }
2628273d9f13SBarry Smith     b->v          = data;
2629273d9f13SBarry Smith     b->user_alloc = PETSC_TRUE;
2630273d9f13SBarry Smith   }
26310450473dSBarry Smith   B->assembled = PETSC_TRUE;
2632273d9f13SBarry Smith   PetscFunctionReturn(0);
2633273d9f13SBarry Smith }
2634273d9f13SBarry Smith 
263565b80a83SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
2636cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat)
26378baccfbdSHong Zhang {
2638d77f618aSHong Zhang   Mat               mat_elemental;
2639d77f618aSHong Zhang   PetscErrorCode    ierr;
26401683a169SBarry Smith   const PetscScalar *array;
26411683a169SBarry Smith   PetscScalar       *v_colwise;
2642d77f618aSHong Zhang   PetscInt          M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols;
2643d77f618aSHong Zhang 
26448baccfbdSHong Zhang   PetscFunctionBegin;
2645d77f618aSHong Zhang   ierr = PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols);CHKERRQ(ierr);
26461683a169SBarry Smith   ierr = MatDenseGetArrayRead(A,&array);CHKERRQ(ierr);
2647d77f618aSHong Zhang   /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */
2648d77f618aSHong Zhang   k = 0;
2649d77f618aSHong Zhang   for (j=0; j<N; j++) {
2650d77f618aSHong Zhang     cols[j] = j;
2651d77f618aSHong Zhang     for (i=0; i<M; i++) {
2652d77f618aSHong Zhang       v_colwise[j*M+i] = array[k++];
2653d77f618aSHong Zhang     }
2654d77f618aSHong Zhang   }
2655d77f618aSHong Zhang   for (i=0; i<M; i++) {
2656d77f618aSHong Zhang     rows[i] = i;
2657d77f618aSHong Zhang   }
26581683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(A,&array);CHKERRQ(ierr);
2659d77f618aSHong Zhang 
2660d77f618aSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental);CHKERRQ(ierr);
2661d77f618aSHong Zhang   ierr = MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
2662d77f618aSHong Zhang   ierr = MatSetType(mat_elemental,MATELEMENTAL);CHKERRQ(ierr);
2663d77f618aSHong Zhang   ierr = MatSetUp(mat_elemental);CHKERRQ(ierr);
2664d77f618aSHong Zhang 
2665d77f618aSHong Zhang   /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */
2666d77f618aSHong Zhang   ierr = MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES);CHKERRQ(ierr);
2667d77f618aSHong Zhang   ierr = MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2668d77f618aSHong Zhang   ierr = MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2669d77f618aSHong Zhang   ierr = PetscFree3(v_colwise,rows,cols);CHKERRQ(ierr);
2670d77f618aSHong Zhang 
2671511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
267228be2f97SBarry Smith     ierr = MatHeaderReplace(A,&mat_elemental);CHKERRQ(ierr);
2673d77f618aSHong Zhang   } else {
2674d77f618aSHong Zhang     *newmat = mat_elemental;
2675d77f618aSHong Zhang   }
26768baccfbdSHong Zhang   PetscFunctionReturn(0);
26778baccfbdSHong Zhang }
267865b80a83SHong Zhang #endif
26798baccfbdSHong Zhang 
26801b807ce4Svictorle /*@C
26811b807ce4Svictorle   MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array
26821b807ce4Svictorle 
26831b807ce4Svictorle   Input parameter:
26841b807ce4Svictorle + A - the matrix
26851b807ce4Svictorle - lda - the leading dimension
26861b807ce4Svictorle 
26871b807ce4Svictorle   Notes:
2688867c911aSBarry Smith   This routine is to be used in conjunction with MatSeqDenseSetPreallocation();
26891b807ce4Svictorle   it asserts that the preallocation has a leading dimension (the LDA parameter
26901b807ce4Svictorle   of Blas and Lapack fame) larger than M, the first dimension of the matrix.
26911b807ce4Svictorle 
26921b807ce4Svictorle   Level: intermediate
26931b807ce4Svictorle 
2694284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize()
2695867c911aSBarry Smith 
26961b807ce4Svictorle @*/
26977087cfbeSBarry Smith PetscErrorCode  MatSeqDenseSetLDA(Mat B,PetscInt lda)
26981b807ce4Svictorle {
26991b807ce4Svictorle   Mat_SeqDense *b = (Mat_SeqDense*)B->data;
270021a2c019SBarry Smith 
27011b807ce4Svictorle   PetscFunctionBegin;
2702e32f2f54SBarry Smith   if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n);
27031b807ce4Svictorle   b->lda       = lda;
270421a2c019SBarry Smith   b->changelda = PETSC_FALSE;
270521a2c019SBarry Smith   b->Mmax      = PetscMax(b->Mmax,lda);
27061b807ce4Svictorle   PetscFunctionReturn(0);
27071b807ce4Svictorle }
27081b807ce4Svictorle 
2709d528f656SJakub Kruzik PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2710d528f656SJakub Kruzik {
2711d528f656SJakub Kruzik   PetscErrorCode ierr;
2712d528f656SJakub Kruzik   PetscMPIInt    size;
2713d528f656SJakub Kruzik 
2714d528f656SJakub Kruzik   PetscFunctionBegin;
2715d528f656SJakub Kruzik   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2716d528f656SJakub Kruzik   if (size == 1) {
2717d528f656SJakub Kruzik     if (scall == MAT_INITIAL_MATRIX) {
2718d528f656SJakub Kruzik       ierr = MatDuplicate(inmat,MAT_COPY_VALUES,outmat);CHKERRQ(ierr);
2719d528f656SJakub Kruzik     } else {
2720d528f656SJakub Kruzik       ierr = MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2721d528f656SJakub Kruzik     }
2722d528f656SJakub Kruzik   } else {
2723d528f656SJakub Kruzik     ierr = MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat);CHKERRQ(ierr);
2724d528f656SJakub Kruzik   }
2725d528f656SJakub Kruzik   PetscFunctionReturn(0);
2726d528f656SJakub Kruzik }
2727d528f656SJakub Kruzik 
27280bad9183SKris Buschelman /*MC
2729fafad747SKris Buschelman    MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices.
27300bad9183SKris Buschelman 
27310bad9183SKris Buschelman    Options Database Keys:
27320bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions()
27330bad9183SKris Buschelman 
27340bad9183SKris Buschelman   Level: beginner
27350bad9183SKris Buschelman 
273689665df3SBarry Smith .seealso: MatCreateSeqDense()
273789665df3SBarry Smith 
27380bad9183SKris Buschelman M*/
27390bad9183SKris Buschelman 
27408cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat B)
2741273d9f13SBarry Smith {
2742273d9f13SBarry Smith   Mat_SeqDense   *b;
2743dfbe8321SBarry Smith   PetscErrorCode ierr;
27447c334f02SBarry Smith   PetscMPIInt    size;
2745273d9f13SBarry Smith 
2746273d9f13SBarry Smith   PetscFunctionBegin;
2747ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr);
2748e32f2f54SBarry Smith   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
274955659b69SBarry Smith 
2750b00a9115SJed Brown   ierr    = PetscNewLog(B,&b);CHKERRQ(ierr);
2751549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
275244cd7ae7SLois Curfman McInnes   B->data = (void*)b;
275318f449edSLois Curfman McInnes 
2754273d9f13SBarry Smith   b->roworiented = PETSC_TRUE;
27554e220ebcSLois Curfman McInnes 
275649a6ff4bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense);CHKERRQ(ierr);
2757bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);CHKERRQ(ierr);
27588572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr);
2759d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense);CHKERRQ(ierr);
2760d3042a70SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense);CHKERRQ(ierr);
27618572280aSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense);CHKERRQ(ierr);
2762715b7558SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr);
2763bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr);
27648baccfbdSHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
27658baccfbdSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental);CHKERRQ(ierr);
27668baccfbdSHong Zhang #endif
2767bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr);
2768bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
2769bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
2770bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2771abc3b08eSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_seqdense_C",MatPtAP_SeqDense_SeqDense);CHKERRQ(ierr);
27724099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijperm_seqdense_C",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
27734099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijperm_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
27744099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijperm_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
27754099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijperm_seqdense_C",MatPtAP_SeqDense_SeqDense);CHKERRQ(ierr);
2776e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijsell_seqdense_C",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
2777e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijsell_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
2778e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijsell_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2779e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijsell_seqdense_C",MatPtAP_SeqDense_SeqDense);CHKERRQ(ierr);
278096e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaijmkl_seqdense_C",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
278196e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaijmkl_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
278296e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaijmkl_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
278396e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaijmkl_seqdense_C",MatPtAP_SeqDense_SeqDense);CHKERRQ(ierr);
278496e6d5c4SRichard Tran Mills 
27853bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaij_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
27863bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
27873bf78175SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
27884099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijperm_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
27894099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijperm_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
27904099cc6bSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijperm_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2791e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijsell_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
2792e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijsell_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
2793e9e4f4a6SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijsell_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
279496e6d5c4SRichard Tran Mills 
279596e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaijmkl_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr);
279696e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaijmkl_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr);
279796e6d5c4SRichard Tran Mills   ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaijmkl_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr);
2798af53bab2SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense);CHKERRQ(ierr);
2799af53bab2SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense);CHKERRQ(ierr);
280017667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr);
28013a40ed3dSBarry Smith   PetscFunctionReturn(0);
2802289bc588SBarry Smith }
280386aefd0dSHong Zhang 
280486aefd0dSHong Zhang /*@C
2805af53bab2SHong Zhang    MatDenseGetColumn - gives access to a column of a dense matrix. This is only the local part of the column. You MUST call MatDenseRestoreColumn() to avoid memory bleeding.
280686aefd0dSHong Zhang 
280786aefd0dSHong Zhang    Not Collective
280886aefd0dSHong Zhang 
280986aefd0dSHong Zhang    Input Parameter:
281086aefd0dSHong Zhang +  mat - a MATSEQDENSE or MATMPIDENSE matrix
281186aefd0dSHong Zhang -  col - column index
281286aefd0dSHong Zhang 
281386aefd0dSHong Zhang    Output Parameter:
281486aefd0dSHong Zhang .  vals - pointer to the data
281586aefd0dSHong Zhang 
281686aefd0dSHong Zhang    Level: intermediate
281786aefd0dSHong Zhang 
281886aefd0dSHong Zhang .seealso: MatDenseRestoreColumn()
281986aefd0dSHong Zhang @*/
282086aefd0dSHong Zhang PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals)
282186aefd0dSHong Zhang {
282286aefd0dSHong Zhang   PetscErrorCode ierr;
282386aefd0dSHong Zhang 
282486aefd0dSHong Zhang   PetscFunctionBegin;
282586aefd0dSHong Zhang   ierr = PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals));CHKERRQ(ierr);
282686aefd0dSHong Zhang   PetscFunctionReturn(0);
282786aefd0dSHong Zhang }
282886aefd0dSHong Zhang 
282986aefd0dSHong Zhang /*@C
283086aefd0dSHong Zhang    MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn().
283186aefd0dSHong Zhang 
283286aefd0dSHong Zhang    Not Collective
283386aefd0dSHong Zhang 
283486aefd0dSHong Zhang    Input Parameter:
283586aefd0dSHong Zhang .  mat - a MATSEQDENSE or MATMPIDENSE matrix
283686aefd0dSHong Zhang 
283786aefd0dSHong Zhang    Output Parameter:
283886aefd0dSHong Zhang .  vals - pointer to the data
283986aefd0dSHong Zhang 
284086aefd0dSHong Zhang    Level: intermediate
284186aefd0dSHong Zhang 
284286aefd0dSHong Zhang .seealso: MatDenseGetColumn()
284386aefd0dSHong Zhang @*/
284486aefd0dSHong Zhang PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals)
284586aefd0dSHong Zhang {
284686aefd0dSHong Zhang   PetscErrorCode ierr;
284786aefd0dSHong Zhang 
284886aefd0dSHong Zhang   PetscFunctionBegin;
284986aefd0dSHong Zhang   ierr = PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals));CHKERRQ(ierr);
285086aefd0dSHong Zhang   PetscFunctionReturn(0);
285186aefd0dSHong Zhang }
2852