1*29bbc08cSBarry Smith /*$Id: dense.c,v 1.188 2000/08/16 18:44:49 balay Exp bsmith $*/ 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 670f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 7b4c862e0SSatish Balay #include "petscblaslapack.h" 8289bc588SBarry Smith 95615d1e5SSatish Balay #undef __FUNC__ 10b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAXPY_SeqDense" 111987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 121987afe7SBarry Smith { 131987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 141987afe7SBarry Smith int N = x->m*x->n,one = 1; 153a40ed3dSBarry Smith 163a40ed3dSBarry Smith PetscFunctionBegin; 171987afe7SBarry Smith BLaxpy_(&N,alpha,x->v,&one,y->v,&one); 181987afe7SBarry Smith PLogFlops(2*N-1); 193a40ed3dSBarry Smith PetscFunctionReturn(0); 201987afe7SBarry Smith } 211987afe7SBarry Smith 225615d1e5SSatish Balay #undef __FUNC__ 23b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetInfo_SeqDense" 248f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 25289bc588SBarry Smith { 264e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 274e220ebcSLois Curfman McInnes int i,N = a->m*a->n,count = 0; 284e220ebcSLois Curfman McInnes Scalar *v = a->v; 293a40ed3dSBarry Smith 303a40ed3dSBarry Smith PetscFunctionBegin; 31289bc588SBarry Smith for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 324e220ebcSLois Curfman McInnes 334e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 344e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 354e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 364e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 374e220ebcSLois Curfman McInnes info->block_size = 1.0; 384e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 394e220ebcSLois Curfman McInnes info->nz_used = (double)count; 404e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 414e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 424e220ebcSLois Curfman McInnes info->mallocs = 0; 434e220ebcSLois Curfman McInnes info->memory = A->mem; 444e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 454e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 464e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 474e220ebcSLois Curfman McInnes 483a40ed3dSBarry Smith PetscFunctionReturn(0); 49289bc588SBarry Smith } 50289bc588SBarry Smith 515615d1e5SSatish Balay #undef __FUNC__ 52b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatScale_SeqDense" 538f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA) 5480cd9d93SLois Curfman McInnes { 5580cd9d93SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)inA->data; 5680cd9d93SLois Curfman McInnes int one = 1,nz; 5780cd9d93SLois Curfman McInnes 583a40ed3dSBarry Smith PetscFunctionBegin; 5980cd9d93SLois Curfman McInnes nz = a->m*a->n; 6080cd9d93SLois Curfman McInnes BLscal_(&nz,alpha,a->v,&one); 6180cd9d93SLois Curfman McInnes PLogFlops(nz); 623a40ed3dSBarry Smith PetscFunctionReturn(0); 6380cd9d93SLois Curfman McInnes } 6480cd9d93SLois Curfman McInnes 65289bc588SBarry Smith /* ---------------------------------------------------------------*/ 666831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 676831982aSBarry Smith rather than put it in the Mat->row slot.*/ 685615d1e5SSatish Balay #undef __FUNC__ 69b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLUFactor_SeqDense" 70e03a110bSBarry Smith int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatLUInfo *minfo) 71289bc588SBarry Smith { 72c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 736abc6512SBarry Smith int info; 7455659b69SBarry Smith 753a40ed3dSBarry Smith PetscFunctionBegin; 76289bc588SBarry Smith if (!mat->pivots) { 7745d48df9SBarry Smith mat->pivots = (int*)PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 78c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 79289bc588SBarry Smith } 80f1af5d2fSBarry Smith A->factor = FACTOR_LU; 817a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 82289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 83*29bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 84*29bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 8555659b69SBarry Smith PLogFlops((2*mat->n*mat->n*mat->n)/3); 863a40ed3dSBarry Smith PetscFunctionReturn(0); 87289bc588SBarry Smith } 886ee01492SSatish Balay 895615d1e5SSatish Balay #undef __FUNC__ 90b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_SeqDense" 915609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9202cad45dSBarry Smith { 9302cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 9402cad45dSBarry Smith int ierr; 9502cad45dSBarry Smith Mat newi; 9602cad45dSBarry Smith 973a40ed3dSBarry Smith PetscFunctionBegin; 9802cad45dSBarry Smith ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi);CHKERRQ(ierr); 9902cad45dSBarry Smith l = (Mat_SeqDense*)newi->data; 1005609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 101549d3d68SSatish Balay ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 10202cad45dSBarry Smith } 10302cad45dSBarry Smith newi->assembled = PETSC_TRUE; 10402cad45dSBarry Smith *newmat = newi; 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 10602cad45dSBarry Smith } 10702cad45dSBarry Smith 1085615d1e5SSatish Balay #undef __FUNC__ 109b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLUFactorSymbolic_SeqDense" 110e03a110bSBarry Smith int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact) 111289bc588SBarry Smith { 1123a40ed3dSBarry Smith int ierr; 1133a40ed3dSBarry Smith 1143a40ed3dSBarry Smith PetscFunctionBegin; 1155609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1163a40ed3dSBarry Smith PetscFunctionReturn(0); 117289bc588SBarry Smith } 1186ee01492SSatish Balay 1195615d1e5SSatish Balay #undef __FUNC__ 120b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLUFactorNumeric_SeqDense" 1218f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 122289bc588SBarry Smith { 12302cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 1243a40ed3dSBarry Smith int ierr; 1253a40ed3dSBarry Smith 1263a40ed3dSBarry Smith PetscFunctionBegin; 12702cad45dSBarry Smith /* copy the numerical values */ 128549d3d68SSatish Balay ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 12902cad45dSBarry Smith (*fact)->factor = 0; 130e03a110bSBarry Smith ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr); 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 132289bc588SBarry Smith } 1336ee01492SSatish Balay 1345615d1e5SSatish Balay #undef __FUNC__ 135b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCholeskyFactorSymbolic_SeqDense" 136329f5518SBarry Smith int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact) 137289bc588SBarry Smith { 1383a40ed3dSBarry Smith int ierr; 1393a40ed3dSBarry Smith 1403a40ed3dSBarry Smith PetscFunctionBegin; 1413a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1423a40ed3dSBarry Smith PetscFunctionReturn(0); 143289bc588SBarry Smith } 1446ee01492SSatish Balay 1455615d1e5SSatish Balay #undef __FUNC__ 146b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCholeskyFactorNumeric_SeqDense" 1478f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 148289bc588SBarry Smith { 1493a40ed3dSBarry Smith int ierr; 1503a40ed3dSBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 1523a40ed3dSBarry Smith ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr); 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 154289bc588SBarry Smith } 1556ee01492SSatish Balay 1565615d1e5SSatish Balay #undef __FUNC__ 1575615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense" 158329f5518SBarry Smith int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f) 159289bc588SBarry Smith { 160c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 161606d414cSSatish Balay int info,ierr; 16255659b69SBarry Smith 1633a40ed3dSBarry Smith PetscFunctionBegin; 164752f5567SLois Curfman McInnes if (mat->pivots) { 165606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 166c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 167752f5567SLois Curfman McInnes mat->pivots = 0; 168752f5567SLois Curfman McInnes } 1697a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 170289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 171*29bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 172c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 17355659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 1743a40ed3dSBarry Smith PetscFunctionReturn(0); 175289bc588SBarry Smith } 176289bc588SBarry Smith 1775615d1e5SSatish Balay #undef __FUNC__ 178b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSolve_SeqDense" 1798f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 180289bc588SBarry Smith { 181c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 182a57ad546SLois Curfman McInnes int one = 1,info,ierr; 1836abc6512SBarry Smith Scalar *x,*y; 18467e560aaSBarry Smith 1853a40ed3dSBarry Smith PetscFunctionBegin; 1867a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 1877a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1887a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 189549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 190c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 19148d91487SBarry Smith LAgetrs_("N",&mat->m,&one,mat->v,&mat->m,mat->pivots,y,&mat->m,&info); 1927a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 19348d91487SBarry Smith LApotrs_("L",&mat->m,&one,mat->v,&mat->m,y,&mat->m,&info); 194289bc588SBarry Smith } 195*29bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 196*29bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 1977a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1987a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 1993e80ec95SBarry Smith PLogFlops(2*mat->n*mat->n - mat->n); 2003a40ed3dSBarry Smith PetscFunctionReturn(0); 201289bc588SBarry Smith } 2026ee01492SSatish Balay 2035615d1e5SSatish Balay #undef __FUNC__ 204b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSolveTranspose_SeqDense" 2057c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 206da3a660dSBarry Smith { 207c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2087a97a34bSBarry Smith int ierr,one = 1,info; 2096abc6512SBarry Smith Scalar *x,*y; 21067e560aaSBarry Smith 2113a40ed3dSBarry Smith PetscFunctionBegin; 2127a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 2137a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2147a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 215549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 216752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 217da3a660dSBarry Smith if (mat->pivots) { 21848d91487SBarry Smith LAgetrs_("T",&mat->m,&one,mat->v,&mat->m,mat->pivots,y,&mat->m,&info); 2197a97a34bSBarry Smith } else { 22048d91487SBarry Smith LApotrs_("L",&mat->m,&one,mat->v,&mat->m,y,&mat->m,&info); 221da3a660dSBarry Smith } 222*29bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 2237a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2247a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 2253e80ec95SBarry Smith PLogFlops(2*mat->n*mat->n - mat->n); 2263a40ed3dSBarry Smith PetscFunctionReturn(0); 227da3a660dSBarry Smith } 2286ee01492SSatish Balay 2295615d1e5SSatish Balay #undef __FUNC__ 230b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSolveAdd_SeqDense" 2318f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 232da3a660dSBarry Smith { 233c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2347a97a34bSBarry Smith int ierr,one = 1,info; 2356abc6512SBarry Smith Scalar *x,*y,sone = 1.0; 236da3a660dSBarry Smith Vec tmp = 0; 23767e560aaSBarry Smith 2383a40ed3dSBarry Smith PetscFunctionBegin; 2397a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2407a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 2417a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 242da3a660dSBarry Smith if (yy == zz) { 24378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 244c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 24578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 246da3a660dSBarry Smith } 247549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 248752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 249da3a660dSBarry Smith if (mat->pivots) { 25048d91487SBarry Smith LAgetrs_("N",&mat->m,&one,mat->v,&mat->m,mat->pivots,y,&mat->m,&info); 251a8c6a408SBarry Smith } else { 25248d91487SBarry Smith LApotrs_("L",&mat->m,&one,mat->v,&mat->m,y,&mat->m,&info); 253da3a660dSBarry Smith } 254*29bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 255a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 256a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 2577a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2587a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 2593e80ec95SBarry Smith PLogFlops(2*mat->n*mat->n); 2603a40ed3dSBarry Smith PetscFunctionReturn(0); 261da3a660dSBarry Smith } 26267e560aaSBarry Smith 2635615d1e5SSatish Balay #undef __FUNC__ 264b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSolveTransposeAdd_SeqDense" 2657c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 266da3a660dSBarry Smith { 267c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2686abc6512SBarry Smith int one = 1,info,ierr; 2696abc6512SBarry Smith Scalar *x,*y,sone = 1.0; 270da3a660dSBarry Smith Vec tmp; 27167e560aaSBarry Smith 2723a40ed3dSBarry Smith PetscFunctionBegin; 2737a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 2747a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2757a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 276da3a660dSBarry Smith if (yy == zz) { 27778b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 278c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 27978b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 280da3a660dSBarry Smith } 281549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 282752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 283da3a660dSBarry Smith if (mat->pivots) { 28448d91487SBarry Smith LAgetrs_("T",&mat->m,&one,mat->v,&mat->m,mat->pivots,y,&mat->m,&info); 2853a40ed3dSBarry Smith } else { 28648d91487SBarry Smith LApotrs_("L",&mat->m,&one,mat->v,&mat->m,y,&mat->m,&info); 287da3a660dSBarry Smith } 288*29bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 28990f02eecSBarry Smith if (tmp) { 29090f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 29190f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 2923a40ed3dSBarry Smith } else { 29390f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 29490f02eecSBarry Smith } 2957a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2967a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 2973e80ec95SBarry Smith PLogFlops(2*mat->n*mat->n); 2983a40ed3dSBarry Smith PetscFunctionReturn(0); 299da3a660dSBarry Smith } 300289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3015615d1e5SSatish Balay #undef __FUNC__ 302b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRelax_SeqDense" 303329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 304329f5518SBarry Smith PetscReal shift,int its,Vec xx) 305289bc588SBarry Smith { 306c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 307289bc588SBarry Smith Scalar *x,*b,*v = mat->v,zero = 0.0,xt; 308bc1b551cSSatish Balay int ierr,m = mat->m,i; 309aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 310bc1b551cSSatish Balay int o = 1; 311bc1b551cSSatish Balay #endif 312289bc588SBarry Smith 3133a40ed3dSBarry Smith PetscFunctionBegin; 314289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3153a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 316bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 317289bc588SBarry Smith } 3187a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3197a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 320289bc588SBarry Smith while (its--) { 321289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 322289bc588SBarry Smith for (i=0; i<m; i++) { 323aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 324f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 325f1747703SBarry Smith not happy about returning a double complex */ 326f1747703SBarry Smith int _i; 327f1747703SBarry Smith Scalar sum = b[i]; 328f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3293f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 330f1747703SBarry Smith } 331f1747703SBarry Smith xt = sum; 332f1747703SBarry Smith #else 333289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 334f1747703SBarry Smith #endif 33555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 336289bc588SBarry Smith } 337289bc588SBarry Smith } 338289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 339289bc588SBarry Smith for (i=m-1; i>=0; i--) { 340aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 341f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 342f1747703SBarry Smith not happy about returning a double complex */ 343f1747703SBarry Smith int _i; 344f1747703SBarry Smith Scalar sum = b[i]; 345f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3463f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 347f1747703SBarry Smith } 348f1747703SBarry Smith xt = sum; 349f1747703SBarry Smith #else 350289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 351f1747703SBarry Smith #endif 35255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 353289bc588SBarry Smith } 354289bc588SBarry Smith } 355289bc588SBarry Smith } 356600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 3577a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3583a40ed3dSBarry Smith PetscFunctionReturn(0); 359289bc588SBarry Smith } 360289bc588SBarry Smith 361289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3625615d1e5SSatish Balay #undef __FUNC__ 363b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTranspose_SeqDense" 3647c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 365289bc588SBarry Smith { 366c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 367289bc588SBarry Smith Scalar *v = mat->v,*x,*y; 3687a97a34bSBarry Smith int ierr,_One=1;Scalar _DOne=1.0,_DZero=0.0; 3693a40ed3dSBarry Smith 3703a40ed3dSBarry Smith PetscFunctionBegin; 3717a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 3727a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3737a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 37448d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 3757a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3767a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 37755659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 379289bc588SBarry Smith } 3806ee01492SSatish Balay 3815615d1e5SSatish Balay #undef __FUNC__ 382b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMult_SeqDense" 38344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 384289bc588SBarry Smith { 385c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 386329f5518SBarry Smith Scalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 387329f5518SBarry Smith int ierr,_One=1; 3883a40ed3dSBarry Smith 3893a40ed3dSBarry Smith PetscFunctionBegin; 3907a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 3917a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3927a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 39348d91487SBarry Smith LAgemv_("N",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 3947a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3957a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 39655659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 3973a40ed3dSBarry Smith PetscFunctionReturn(0); 398289bc588SBarry Smith } 3996ee01492SSatish Balay 4005615d1e5SSatish Balay #undef __FUNC__ 401b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultAdd_SeqDense" 40244cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 403289bc588SBarry Smith { 404c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 405329f5518SBarry Smith Scalar *v = mat->v,*x,*y,_DOne=1.0; 406329f5518SBarry Smith int ierr,_One=1; 4073a40ed3dSBarry Smith 4083a40ed3dSBarry Smith PetscFunctionBegin; 4097a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 410600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4117a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4127a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 41348d91487SBarry Smith LAgemv_("N",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 4147a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4157a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 41655659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 4173a40ed3dSBarry Smith PetscFunctionReturn(0); 418289bc588SBarry Smith } 4196ee01492SSatish Balay 4205615d1e5SSatish Balay #undef __FUNC__ 421b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMultTransposeAdd_SeqDense" 4227c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 423289bc588SBarry Smith { 424c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 425600d6d0bSBarry Smith Scalar *v = mat->v,*x,*y; 4267a97a34bSBarry Smith int ierr,_One=1; 4276abc6512SBarry Smith Scalar _DOne=1.0; 4283a40ed3dSBarry Smith 4293a40ed3dSBarry Smith PetscFunctionBegin; 4307a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 431600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4327a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4337a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 43448d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 4357a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4367a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 43755659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 4383a40ed3dSBarry Smith PetscFunctionReturn(0); 439289bc588SBarry Smith } 440289bc588SBarry Smith 441289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4425615d1e5SSatish Balay #undef __FUNC__ 443b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_SeqDense" 4448f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 445289bc588SBarry Smith { 446c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 447289bc588SBarry Smith Scalar *v; 448289bc588SBarry Smith int i; 44967e560aaSBarry Smith 4503a40ed3dSBarry Smith PetscFunctionBegin; 451289bc588SBarry Smith *ncols = mat->n; 452289bc588SBarry Smith if (cols) { 45345d48df9SBarry Smith *cols = (int*)PetscMalloc((mat->n+1)*sizeof(int));CHKPTRQ(*cols); 45473c3ce57SBarry Smith for (i=0; i<mat->n; i++) (*cols)[i] = i; 455289bc588SBarry Smith } 456289bc588SBarry Smith if (vals) { 45745d48df9SBarry Smith *vals = (Scalar*)PetscMalloc((mat->n+1)*sizeof(Scalar));CHKPTRQ(*vals); 458289bc588SBarry Smith v = mat->v + row; 45973c3ce57SBarry Smith for (i=0; i<mat->n; i++) {(*vals)[i] = *v; v += mat->m;} 460289bc588SBarry Smith } 4613a40ed3dSBarry Smith PetscFunctionReturn(0); 462289bc588SBarry Smith } 4636ee01492SSatish Balay 4645615d1e5SSatish Balay #undef __FUNC__ 465b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_SeqDense" 4668f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 467289bc588SBarry Smith { 468606d414cSSatish Balay int ierr; 469606d414cSSatish Balay PetscFunctionBegin; 470606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 471606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 4723a40ed3dSBarry Smith PetscFunctionReturn(0); 473289bc588SBarry Smith } 474289bc588SBarry Smith /* ----------------------------------------------------------------*/ 4755615d1e5SSatish Balay #undef __FUNC__ 476b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValues_SeqDense" 4778f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 478e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 479289bc588SBarry Smith { 480c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 481289bc588SBarry Smith int i,j; 482d6dfbf8fSBarry Smith 4833a40ed3dSBarry Smith PetscFunctionBegin; 484289bc588SBarry Smith if (!mat->roworiented) { 485dbb450caSBarry Smith if (addv == INSERT_VALUES) { 486289bc588SBarry Smith for (j=0; j<n; j++) { 4875ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 488aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 489*29bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 49058804f6dSBarry Smith #endif 491289bc588SBarry Smith for (i=0; i<m; i++) { 4925ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 493aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 494*29bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 49558804f6dSBarry Smith #endif 496289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 497289bc588SBarry Smith } 498289bc588SBarry Smith } 4993a40ed3dSBarry Smith } else { 500289bc588SBarry Smith for (j=0; j<n; j++) { 5015ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 502aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 503*29bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 50458804f6dSBarry Smith #endif 505289bc588SBarry Smith for (i=0; i<m; i++) { 5065ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 507aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 508*29bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 50958804f6dSBarry Smith #endif 510289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 511289bc588SBarry Smith } 512289bc588SBarry Smith } 513289bc588SBarry Smith } 5143a40ed3dSBarry Smith } else { 515dbb450caSBarry Smith if (addv == INSERT_VALUES) { 516e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 5175ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 518aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 519*29bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 52058804f6dSBarry Smith #endif 521e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 5225ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 523aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 524*29bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 52558804f6dSBarry Smith #endif 526e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 527e8d4e0b9SBarry Smith } 528e8d4e0b9SBarry Smith } 5293a40ed3dSBarry Smith } else { 530289bc588SBarry Smith for (i=0; i<m; i++) { 5315ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 532aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 533*29bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 53458804f6dSBarry Smith #endif 535289bc588SBarry Smith for (j=0; j<n; j++) { 5365ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 537aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 538*29bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 53958804f6dSBarry Smith #endif 540289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 541289bc588SBarry Smith } 542289bc588SBarry Smith } 543289bc588SBarry Smith } 544e8d4e0b9SBarry Smith } 5453a40ed3dSBarry Smith PetscFunctionReturn(0); 546289bc588SBarry Smith } 547e8d4e0b9SBarry Smith 5485615d1e5SSatish Balay #undef __FUNC__ 549b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetValues_SeqDense" 5508f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 551ae80bb75SLois Curfman McInnes { 552ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 553ae80bb75SLois Curfman McInnes int i,j; 554ae80bb75SLois Curfman McInnes Scalar *vpt = v; 555ae80bb75SLois Curfman McInnes 5563a40ed3dSBarry Smith PetscFunctionBegin; 557ae80bb75SLois Curfman McInnes /* row-oriented output */ 558ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 559ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 560ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 561ae80bb75SLois Curfman McInnes } 562ae80bb75SLois Curfman McInnes } 5633a40ed3dSBarry Smith PetscFunctionReturn(0); 564ae80bb75SLois Curfman McInnes } 565ae80bb75SLois Curfman McInnes 566289bc588SBarry Smith /* -----------------------------------------------------------------*/ 567289bc588SBarry Smith 568e090d566SSatish Balay #include "petscsys.h" 56988e32edaSLois Curfman McInnes 5705615d1e5SSatish Balay #undef __FUNC__ 571b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_SeqDense" 57219bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 57388e32edaSLois Curfman McInnes { 57488e32edaSLois Curfman McInnes Mat_SeqDense *a; 57588e32edaSLois Curfman McInnes Mat B; 57688e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 57788e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 5784a41a4d3SSatish Balay Scalar *vals,*svals,*v,*w; 57919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 58088e32edaSLois Curfman McInnes 5813a40ed3dSBarry Smith PetscFunctionBegin; 582d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 583*29bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 58490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 5850752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 586*29bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 58788e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 58888e32edaSLois Curfman McInnes 58990ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 59090ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 59190ace30eSBarry Smith B = *A; 59290ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 5934a41a4d3SSatish Balay v = a->v; 5944a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 5954a41a4d3SSatish Balay from row major to column major */ 5964a41a4d3SSatish Balay w = (Scalar*)PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w); 59790ace30eSBarry Smith /* read in nonzero values */ 5984a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 5994a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6004a41a4d3SSatish Balay for (j=0; j<N; j++) { 6014a41a4d3SSatish Balay for (i=0; i<M; i++) { 6024a41a4d3SSatish Balay *v++ =w[i*N+j]; 6034a41a4d3SSatish Balay } 6044a41a4d3SSatish Balay } 605606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6066d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6076d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60890ace30eSBarry Smith } else { 60988e32edaSLois Curfman McInnes /* read row lengths */ 61045d48df9SBarry Smith rowlengths = (int*)PetscMalloc((M+1)*sizeof(int));CHKPTRQ(rowlengths); 6110752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 61288e32edaSLois Curfman McInnes 61388e32edaSLois Curfman McInnes /* create our matrix */ 614b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 61588e32edaSLois Curfman McInnes B = *A; 61688e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 61788e32edaSLois Curfman McInnes v = a->v; 61888e32edaSLois Curfman McInnes 61988e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 62045d48df9SBarry Smith cols = scols = (int*)PetscMalloc((nz+1)*sizeof(int));CHKPTRQ(cols); 6210752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 62245d48df9SBarry Smith vals = svals = (Scalar*)PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(vals); 6230752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 62488e32edaSLois Curfman McInnes 62588e32edaSLois Curfman McInnes /* insert into matrix */ 62688e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 62788e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 62888e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 62988e32edaSLois Curfman McInnes } 630606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 631606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 632606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 63388e32edaSLois Curfman McInnes 6346d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6356d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63690ace30eSBarry Smith } 6373a40ed3dSBarry Smith PetscFunctionReturn(0); 63888e32edaSLois Curfman McInnes } 63988e32edaSLois Curfman McInnes 640e090d566SSatish Balay #include "petscsys.h" 641932b0c3eSLois Curfman McInnes 6425615d1e5SSatish Balay #undef __FUNC__ 643b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_ASCII" 644932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 645289bc588SBarry Smith { 646932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 647932b0c3eSLois Curfman McInnes int ierr,i,j,format; 648932b0c3eSLois Curfman McInnes char *outputname; 649932b0c3eSLois Curfman McInnes Scalar *v; 650932b0c3eSLois Curfman McInnes 6513a40ed3dSBarry Smith PetscFunctionBegin; 65277ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 653888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 654639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 6553a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 6566831982aSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 6576831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 65880cd9d93SLois Curfman McInnes for (i=0; i<a->m; i++) { 65944cd7ae7SLois Curfman McInnes v = a->v + i; 6606831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 66180cd9d93SLois Curfman McInnes for (j=0; j<a->n; j++) { 662aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 663329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 664329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 665329f5518SBarry Smith } else if (PetscRealPart(*v)) { 666329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));CHKERRQ(ierr); 6676831982aSBarry Smith } 66880cd9d93SLois Curfman McInnes #else 6696831982aSBarry Smith if (*v) { 6706831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr); 6716831982aSBarry Smith } 67280cd9d93SLois Curfman McInnes #endif 673d64ed03dSBarry Smith v += a->m; 67480cd9d93SLois Curfman McInnes } 6756831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 67680cd9d93SLois Curfman McInnes } 6776831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 6783a40ed3dSBarry Smith } else { 6796831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 680aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 68147989497SBarry Smith int allreal = 1; 68247989497SBarry Smith /* determine if matrix has all real values */ 68347989497SBarry Smith v = a->v; 68447989497SBarry Smith for (i=0; i<a->m*a->n; i++) { 685329f5518SBarry Smith if (PetscImaginaryPart(v[i])) { allreal = 0; break ;} 68647989497SBarry Smith } 68747989497SBarry Smith #endif 688932b0c3eSLois Curfman McInnes for (i=0; i<a->m; i++) { 689932b0c3eSLois Curfman McInnes v = a->v + i; 690932b0c3eSLois Curfman McInnes for (j=0; j<a->n; j++) { 691aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 69247989497SBarry Smith if (allreal) { 693329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += a->m; 69447989497SBarry Smith } else { 695329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += a->m; 69647989497SBarry Smith } 697289bc588SBarry Smith #else 6986831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += a->m; 699289bc588SBarry Smith #endif 700289bc588SBarry Smith } 7016831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 702289bc588SBarry Smith } 7036831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 704da3a660dSBarry Smith } 7056831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 7063a40ed3dSBarry Smith PetscFunctionReturn(0); 707289bc588SBarry Smith } 708289bc588SBarry Smith 7095615d1e5SSatish Balay #undef __FUNC__ 710b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Binary" 711932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 712932b0c3eSLois Curfman McInnes { 713932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 714932b0c3eSLois Curfman McInnes int ict,j,n = a->n,m = a->m,i,fd,*col_lens,ierr,nz = m*n; 71590ace30eSBarry Smith int format; 71690ace30eSBarry Smith Scalar *v,*anonz,*vals; 717932b0c3eSLois Curfman McInnes 7183a40ed3dSBarry Smith PetscFunctionBegin; 71990ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 72090ace30eSBarry Smith 72190ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 72202cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 72390ace30eSBarry Smith /* store the matrix as a dense matrix */ 72490ace30eSBarry Smith col_lens = (int*)PetscMalloc(4*sizeof(int));CHKPTRQ(col_lens); 72590ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 72690ace30eSBarry Smith col_lens[1] = m; 72790ace30eSBarry Smith col_lens[2] = n; 72890ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 7290752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 730606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 73190ace30eSBarry Smith 73290ace30eSBarry Smith /* write out matrix, by rows */ 73345d48df9SBarry Smith vals = (Scalar*)PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(vals); 73490ace30eSBarry Smith v = a->v; 73590ace30eSBarry Smith for (i=0; i<m; i++) { 73690ace30eSBarry Smith for (j=0; j<n; j++) { 73790ace30eSBarry Smith vals[i + j*m] = *v++; 73890ace30eSBarry Smith } 73990ace30eSBarry Smith } 7400752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 741606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 74290ace30eSBarry Smith } else { 7430452661fSBarry Smith col_lens = (int*)PetscMalloc((4+nz)*sizeof(int));CHKPTRQ(col_lens); 744932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 745932b0c3eSLois Curfman McInnes col_lens[1] = m; 746932b0c3eSLois Curfman McInnes col_lens[2] = n; 747932b0c3eSLois Curfman McInnes col_lens[3] = nz; 748932b0c3eSLois Curfman McInnes 749932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 750932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 7510752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 752932b0c3eSLois Curfman McInnes 753932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 754932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 755932b0c3eSLois Curfman McInnes ict = 0; 756932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 757932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 758932b0c3eSLois Curfman McInnes } 7590752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 760606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 761932b0c3eSLois Curfman McInnes 762932b0c3eSLois Curfman McInnes /* store nonzero values */ 76345d48df9SBarry Smith anonz = (Scalar*)PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(anonz); 764932b0c3eSLois Curfman McInnes ict = 0; 765932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 766932b0c3eSLois Curfman McInnes v = a->v + i; 767932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 768932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 769932b0c3eSLois Curfman McInnes } 770932b0c3eSLois Curfman McInnes } 7710752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 772606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 77390ace30eSBarry Smith } 7743a40ed3dSBarry Smith PetscFunctionReturn(0); 775932b0c3eSLois Curfman McInnes } 776932b0c3eSLois Curfman McInnes 7775615d1e5SSatish Balay #undef __FUNC__ 778b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Draw_Zoom" 779f1af5d2fSBarry Smith int MatView_SeqDense_Draw_Zoom(Draw draw,void *Aa) 780f1af5d2fSBarry Smith { 781f1af5d2fSBarry Smith Mat A = (Mat) Aa; 782f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 783f1af5d2fSBarry Smith int m = a->m,n = a->n,format,color,i,j,ierr; 784f1af5d2fSBarry Smith Scalar *v = a->v; 785f1af5d2fSBarry Smith Viewer viewer; 786f1af5d2fSBarry Smith Draw popup; 787329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 788f1af5d2fSBarry Smith 789f1af5d2fSBarry Smith PetscFunctionBegin; 790f1af5d2fSBarry Smith 791f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 792f1af5d2fSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 793f1af5d2fSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 794f1af5d2fSBarry Smith 795f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 796f1af5d2fSBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 797f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 798f1af5d2fSBarry Smith color = DRAW_BLUE; 799f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 800f1af5d2fSBarry Smith x_l = j; 801f1af5d2fSBarry Smith x_r = x_l + 1.0; 802f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 803f1af5d2fSBarry Smith y_l = m - i - 1.0; 804f1af5d2fSBarry Smith y_r = y_l + 1.0; 805f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 806329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 807f1af5d2fSBarry Smith color = DRAW_RED; 808329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 809f1af5d2fSBarry Smith color = DRAW_BLUE; 810f1af5d2fSBarry Smith } else { 811f1af5d2fSBarry Smith continue; 812f1af5d2fSBarry Smith } 813f1af5d2fSBarry Smith #else 814f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 815f1af5d2fSBarry Smith color = DRAW_RED; 816f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 817f1af5d2fSBarry Smith color = DRAW_BLUE; 818f1af5d2fSBarry Smith } else { 819f1af5d2fSBarry Smith continue; 820f1af5d2fSBarry Smith } 821f1af5d2fSBarry Smith #endif 822f1af5d2fSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 823f1af5d2fSBarry Smith } 824f1af5d2fSBarry Smith } 825f1af5d2fSBarry Smith } else { 826f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 827f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 828f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 829f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 830f1af5d2fSBarry Smith } 831f1af5d2fSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 832f1af5d2fSBarry Smith ierr = DrawGetPopup(draw,&popup);CHKERRQ(ierr); 8337c922b88SBarry Smith if (popup) {ierr = DrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 834f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 835f1af5d2fSBarry Smith x_l = j; 836f1af5d2fSBarry Smith x_r = x_l + 1.0; 837f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 838f1af5d2fSBarry Smith y_l = m - i - 1.0; 839f1af5d2fSBarry Smith y_r = y_l + 1.0; 840f1af5d2fSBarry Smith color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 841f1af5d2fSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 842f1af5d2fSBarry Smith } 843f1af5d2fSBarry Smith } 844f1af5d2fSBarry Smith } 845f1af5d2fSBarry Smith PetscFunctionReturn(0); 846f1af5d2fSBarry Smith } 847f1af5d2fSBarry Smith 848f1af5d2fSBarry Smith #undef __FUNC__ 849b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Draw" 850f1af5d2fSBarry Smith int MatView_SeqDense_Draw(Mat A,Viewer viewer) 851f1af5d2fSBarry Smith { 852f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 853f1af5d2fSBarry Smith Draw draw; 854f1af5d2fSBarry Smith PetscTruth isnull; 855329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 856f1af5d2fSBarry Smith int ierr; 857f1af5d2fSBarry Smith 858f1af5d2fSBarry Smith PetscFunctionBegin; 859f1af5d2fSBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 860f1af5d2fSBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); 861f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 862f1af5d2fSBarry Smith 863f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 864f1af5d2fSBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 865f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 866f1af5d2fSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 867f1af5d2fSBarry Smith ierr = DrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 868f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 869f1af5d2fSBarry Smith PetscFunctionReturn(0); 870f1af5d2fSBarry Smith } 871f1af5d2fSBarry Smith 872f1af5d2fSBarry Smith #undef __FUNC__ 873b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense" 874c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer) 875932b0c3eSLois Curfman McInnes { 876932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 877bcd2baecSBarry Smith int ierr; 878f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 879932b0c3eSLois Curfman McInnes 8803a40ed3dSBarry Smith PetscFunctionBegin; 8816831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 8826831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8836831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 884f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 8850f5bd95cSBarry Smith 8860f5bd95cSBarry Smith if (issocket) { 887e03a110bSBarry Smith ierr = ViewerSocketPutScalar(viewer,a->m,a->n,a->v);CHKERRQ(ierr); 8880f5bd95cSBarry Smith } else if (isascii) { 8893a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 8900f5bd95cSBarry Smith } else if (isbinary) { 8913a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 892f1af5d2fSBarry Smith } else if (isdraw) { 893f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 8945cd90555SBarry Smith } else { 895*29bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 896932b0c3eSLois Curfman McInnes } 8973a40ed3dSBarry Smith PetscFunctionReturn(0); 898932b0c3eSLois Curfman McInnes } 899289bc588SBarry Smith 9005615d1e5SSatish Balay #undef __FUNC__ 901b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_SeqDense" 902c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 903289bc588SBarry Smith { 904ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 90590f02eecSBarry Smith int ierr; 90690f02eecSBarry Smith 9073a40ed3dSBarry Smith PetscFunctionBegin; 90894d884c6SBarry Smith 90994d884c6SBarry Smith if (mat->mapping) { 91094d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 91194d884c6SBarry Smith } 91294d884c6SBarry Smith if (mat->bmapping) { 91394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 91494d884c6SBarry Smith } 915aa482453SBarry Smith #if defined(PETSC_USE_LOG) 916c6e7dd08SBarry Smith PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n); 917a5a9c739SBarry Smith #endif 918606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 919606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 920606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 921a5ae1ecdSBarry Smith if (mat->rmap) { 922a5ae1ecdSBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 923a5ae1ecdSBarry Smith } 924a5ae1ecdSBarry Smith if (mat->cmap) { 925a5ae1ecdSBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 92690f02eecSBarry Smith } 927a5a9c739SBarry Smith PLogObjectDestroy(mat); 9280452661fSBarry Smith PetscHeaderDestroy(mat); 9293a40ed3dSBarry Smith PetscFunctionReturn(0); 930289bc588SBarry Smith } 931289bc588SBarry Smith 9325615d1e5SSatish Balay #undef __FUNC__ 933b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_SeqDense" 9348f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 935289bc588SBarry Smith { 936c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 937549d3d68SSatish Balay int k,j,m,n,ierr; 938d3e5ee88SLois Curfman McInnes Scalar *v,tmp; 93948b35521SBarry Smith 9403a40ed3dSBarry Smith PetscFunctionBegin; 941d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 9427c922b88SBarry Smith if (!matout) { /* in place transpose */ 943d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 944d64ed03dSBarry Smith Scalar *w = (Scalar*)PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 945d64ed03dSBarry Smith for (j=0; j<m; j++) { 946d64ed03dSBarry Smith for (k=0; k<n; k++) { 947d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 948d64ed03dSBarry Smith } 949d64ed03dSBarry Smith } 950549d3d68SSatish Balay ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr); 951606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 952d64ed03dSBarry Smith } else { 953d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 954289bc588SBarry Smith for (k=0; k<j; k++) { 955d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 956d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 957d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 958289bc588SBarry Smith } 959289bc588SBarry Smith } 960d64ed03dSBarry Smith } 9613a40ed3dSBarry Smith } else { /* out-of-place transpose */ 962d3e5ee88SLois Curfman McInnes Mat tmat; 963ec8511deSBarry Smith Mat_SeqDense *tmatd; 964d3e5ee88SLois Curfman McInnes Scalar *v2; 965b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 966ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 9670de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 968d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 9690de55854SLois Curfman McInnes for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m]; 970d3e5ee88SLois Curfman McInnes } 9716d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9726d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 973d3e5ee88SLois Curfman McInnes *matout = tmat; 97448b35521SBarry Smith } 9753a40ed3dSBarry Smith PetscFunctionReturn(0); 976289bc588SBarry Smith } 977289bc588SBarry Smith 9785615d1e5SSatish Balay #undef __FUNC__ 979b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatEqual_SeqDense" 9808f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 981289bc588SBarry Smith { 982c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 983c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 984289bc588SBarry Smith int i; 985289bc588SBarry Smith Scalar *v1 = mat1->v,*v2 = mat2->v; 9869ea5d5aeSSatish Balay 9873a40ed3dSBarry Smith PetscFunctionBegin; 988*29bbc08cSBarry Smith if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 9893a40ed3dSBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 9903a40ed3dSBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 991289bc588SBarry Smith for (i=0; i<mat1->m*mat1->n; i++) { 9923a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 993289bc588SBarry Smith v1++; v2++; 994289bc588SBarry Smith } 99577c4ece6SBarry Smith *flg = PETSC_TRUE; 9963a40ed3dSBarry Smith PetscFunctionReturn(0); 997289bc588SBarry Smith } 998289bc588SBarry Smith 9995615d1e5SSatish Balay #undef __FUNC__ 1000b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetDiagonal_SeqDense" 10018f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1002289bc588SBarry Smith { 1003c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10047a97a34bSBarry Smith int ierr,i,n,len; 100544cd7ae7SLois Curfman McInnes Scalar *x,zero = 0.0; 100644cd7ae7SLois Curfman McInnes 10073a40ed3dSBarry Smith PetscFunctionBegin; 10087a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10097a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1010600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 101144cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 1012*29bbc08cSBarry Smith if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 101344cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 1014289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 1015289bc588SBarry Smith } 10167a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10173a40ed3dSBarry Smith PetscFunctionReturn(0); 1018289bc588SBarry Smith } 1019289bc588SBarry Smith 10205615d1e5SSatish Balay #undef __FUNC__ 1021b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDiagonalScale_SeqDense" 10228f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1023289bc588SBarry Smith { 1024c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1025da3a660dSBarry Smith Scalar *l,*r,x,*v; 10267a97a34bSBarry Smith int ierr,i,j,m = mat->m,n = mat->n; 102755659b69SBarry Smith 10283a40ed3dSBarry Smith PetscFunctionBegin; 102928988994SBarry Smith if (ll) { 10307a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1031600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1032*29bbc08cSBarry Smith if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1033da3a660dSBarry Smith for (i=0; i<m; i++) { 1034da3a660dSBarry Smith x = l[i]; 1035da3a660dSBarry Smith v = mat->v + i; 1036da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1037da3a660dSBarry Smith } 10387a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 103944cd7ae7SLois Curfman McInnes PLogFlops(n*m); 1040da3a660dSBarry Smith } 104128988994SBarry Smith if (rr) { 10427a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1043600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1044*29bbc08cSBarry Smith if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1045da3a660dSBarry Smith for (i=0; i<n; i++) { 1046da3a660dSBarry Smith x = r[i]; 1047da3a660dSBarry Smith v = mat->v + i*m; 1048da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1049da3a660dSBarry Smith } 10507a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 105144cd7ae7SLois Curfman McInnes PLogFlops(n*m); 1052da3a660dSBarry Smith } 10533a40ed3dSBarry Smith PetscFunctionReturn(0); 1054289bc588SBarry Smith } 1055289bc588SBarry Smith 10565615d1e5SSatish Balay #undef __FUNC__ 1057b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatNorm_SeqDense" 1058329f5518SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm) 1059289bc588SBarry Smith { 1060c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1061289bc588SBarry Smith Scalar *v = mat->v; 1062329f5518SBarry Smith PetscReal sum = 0.0; 1063289bc588SBarry Smith int i,j; 106455659b69SBarry Smith 10653a40ed3dSBarry Smith PetscFunctionBegin; 1066289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1067289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++) { 1068aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1069329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1070289bc588SBarry Smith #else 1071289bc588SBarry Smith sum += (*v)*(*v); v++; 1072289bc588SBarry Smith #endif 1073289bc588SBarry Smith } 1074289bc588SBarry Smith *norm = sqrt(sum); 107555659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 10763a40ed3dSBarry Smith } else if (type == NORM_1) { 1077289bc588SBarry Smith *norm = 0.0; 1078289bc588SBarry Smith for (j=0; j<mat->n; j++) { 1079289bc588SBarry Smith sum = 0.0; 1080289bc588SBarry Smith for (i=0; i<mat->m; i++) { 108133a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1082289bc588SBarry Smith } 1083289bc588SBarry Smith if (sum > *norm) *norm = sum; 1084289bc588SBarry Smith } 108555659b69SBarry Smith PLogFlops(mat->n*mat->m); 10863a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1087289bc588SBarry Smith *norm = 0.0; 1088289bc588SBarry Smith for (j=0; j<mat->m; j++) { 1089289bc588SBarry Smith v = mat->v + j; 1090289bc588SBarry Smith sum = 0.0; 1091289bc588SBarry Smith for (i=0; i<mat->n; i++) { 1092cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 1093289bc588SBarry Smith } 1094289bc588SBarry Smith if (sum > *norm) *norm = sum; 1095289bc588SBarry Smith } 109655659b69SBarry Smith PLogFlops(mat->n*mat->m); 10973a40ed3dSBarry Smith } else { 1098*29bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1099289bc588SBarry Smith } 11003a40ed3dSBarry Smith PetscFunctionReturn(0); 1101289bc588SBarry Smith } 1102289bc588SBarry Smith 11035615d1e5SSatish Balay #undef __FUNC__ 1104b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_SeqDense" 11058f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1106289bc588SBarry Smith { 1107c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 110867e560aaSBarry Smith 11093a40ed3dSBarry Smith PetscFunctionBegin; 11106d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 11116d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 11126d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 1113219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 11146d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 1115219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 11166d4a8577SBarry Smith op == MAT_SYMMETRIC || 11176d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 11186d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 11196d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 11204787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 11216d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 112290f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 1123b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 1124c38d4ed2SBarry Smith op == MAT_USE_HASH_TABLE) { 1125981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1126c38d4ed2SBarry Smith } else { 1127*29bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 11283a40ed3dSBarry Smith } 11293a40ed3dSBarry Smith PetscFunctionReturn(0); 1130289bc588SBarry Smith } 1131289bc588SBarry Smith 11325615d1e5SSatish Balay #undef __FUNC__ 1133b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroEntries_SeqDense" 11348f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 11356f0a148fSBarry Smith { 1136ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1137549d3d68SSatish Balay int ierr; 11383a40ed3dSBarry Smith 11393a40ed3dSBarry Smith PetscFunctionBegin; 1140549d3d68SSatish Balay ierr = PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));CHKERRQ(ierr); 11413a40ed3dSBarry Smith PetscFunctionReturn(0); 11426f0a148fSBarry Smith } 11436f0a148fSBarry Smith 11445615d1e5SSatish Balay #undef __FUNC__ 1145b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_SeqDense" 11468f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 11474e220ebcSLois Curfman McInnes { 11483a40ed3dSBarry Smith PetscFunctionBegin; 11494e220ebcSLois Curfman McInnes *bs = 1; 11503a40ed3dSBarry Smith PetscFunctionReturn(0); 11514e220ebcSLois Curfman McInnes } 11524e220ebcSLois Curfman McInnes 11535615d1e5SSatish Balay #undef __FUNC__ 1154b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqDense" 11558f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 11566f0a148fSBarry Smith { 1157ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 11586abc6512SBarry Smith int n = l->n,i,j,ierr,N,*rows; 11596f0a148fSBarry Smith Scalar *slot; 116055659b69SBarry Smith 11613a40ed3dSBarry Smith PetscFunctionBegin; 1162e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 116378b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 11646f0a148fSBarry Smith for (i=0; i<N; i++) { 11656f0a148fSBarry Smith slot = l->v + rows[i]; 11666f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 11676f0a148fSBarry Smith } 11686f0a148fSBarry Smith if (diag) { 11696f0a148fSBarry Smith for (i=0; i<N; i++) { 11706f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1171260925b8SBarry Smith *slot = *diag; 11726f0a148fSBarry Smith } 11736f0a148fSBarry Smith } 1174260925b8SBarry Smith ISRestoreIndices(is,&rows); 11753a40ed3dSBarry Smith PetscFunctionReturn(0); 11766f0a148fSBarry Smith } 1177557bce09SLois Curfman McInnes 11785615d1e5SSatish Balay #undef __FUNC__ 1179b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_SeqDense" 11808f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n) 1181557bce09SLois Curfman McInnes { 1182c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11833a40ed3dSBarry Smith 11843a40ed3dSBarry Smith PetscFunctionBegin; 1185f1af5d2fSBarry Smith if (m) *m = mat->m; 1186f1af5d2fSBarry Smith if (n) *n = mat->n; 11873a40ed3dSBarry Smith PetscFunctionReturn(0); 1188557bce09SLois Curfman McInnes } 1189557bce09SLois Curfman McInnes 11905615d1e5SSatish Balay #undef __FUNC__ 1191b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_SeqDense" 11928f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1193d3e5ee88SLois Curfman McInnes { 1194c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11953a40ed3dSBarry Smith 11963a40ed3dSBarry Smith PetscFunctionBegin; 11974c49b128SBarry Smith if (m) *m = 0; 11984c49b128SBarry Smith if (n) *n = mat->m; 11993a40ed3dSBarry Smith PetscFunctionReturn(0); 1200d3e5ee88SLois Curfman McInnes } 1201d3e5ee88SLois Curfman McInnes 12025615d1e5SSatish Balay #undef __FUNC__ 1203b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetArray_SeqDense" 12048f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 120564e87e97SBarry Smith { 1206c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12073a40ed3dSBarry Smith 12083a40ed3dSBarry Smith PetscFunctionBegin; 120964e87e97SBarry Smith *array = mat->v; 12103a40ed3dSBarry Smith PetscFunctionReturn(0); 121164e87e97SBarry Smith } 12120754003eSLois Curfman McInnes 12135615d1e5SSatish Balay #undef __FUNC__ 1214b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreArray_SeqDense" 12158f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1216ff14e315SSatish Balay { 12173a40ed3dSBarry Smith PetscFunctionBegin; 121809b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12193a40ed3dSBarry Smith PetscFunctionReturn(0); 1220ff14e315SSatish Balay } 12210754003eSLois Curfman McInnes 12225615d1e5SSatish Balay #undef __FUNC__ 1223b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrix_SeqDense" 12247b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12250754003eSLois Curfman McInnes { 1226c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1227182d2002SSatish Balay int i,j,ierr,m = mat->m,*irow,*icol,nrows,ncols; 1228182d2002SSatish Balay Scalar *av,*bv,*v = mat->v; 12290754003eSLois Curfman McInnes Mat newmat; 12300754003eSLois Curfman McInnes 12313a40ed3dSBarry Smith PetscFunctionBegin; 123278b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 123378b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1234e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1235e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 12360754003eSLois Curfman McInnes 1237182d2002SSatish Balay /* Check submatrixcall */ 1238182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1239182d2002SSatish Balay int n_cols,n_rows; 1240182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1241*29bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1242182d2002SSatish Balay newmat = *B; 1243182d2002SSatish Balay } else { 12440754003eSLois Curfman McInnes /* Create and fill new matrix */ 1245b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1246182d2002SSatish Balay } 1247182d2002SSatish Balay 1248182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1249182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1250182d2002SSatish Balay 1251182d2002SSatish Balay for (i=0; i<ncols; i++) { 1252182d2002SSatish Balay av = v + m*icol[i]; 1253182d2002SSatish Balay for (j=0; j<nrows; j++) { 1254182d2002SSatish Balay *bv++ = av[irow[j]]; 12550754003eSLois Curfman McInnes } 12560754003eSLois Curfman McInnes } 1257182d2002SSatish Balay 1258182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 12596d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12606d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12610754003eSLois Curfman McInnes 12620754003eSLois Curfman McInnes /* Free work space */ 126378b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 126478b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1265182d2002SSatish Balay *B = newmat; 12663a40ed3dSBarry Smith PetscFunctionReturn(0); 12670754003eSLois Curfman McInnes } 12680754003eSLois Curfman McInnes 12695615d1e5SSatish Balay #undef __FUNC__ 1270b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_SeqDense" 12717b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1272905e6a2fSBarry Smith { 1273905e6a2fSBarry Smith int ierr,i; 1274905e6a2fSBarry Smith 12753a40ed3dSBarry Smith PetscFunctionBegin; 1276905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1277905e6a2fSBarry Smith *B = (Mat*)PetscMalloc((n+1)*sizeof(Mat));CHKPTRQ(*B); 1278905e6a2fSBarry Smith } 1279905e6a2fSBarry Smith 1280905e6a2fSBarry Smith for (i=0; i<n; i++) { 12816a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1282905e6a2fSBarry Smith } 12833a40ed3dSBarry Smith PetscFunctionReturn(0); 1284905e6a2fSBarry Smith } 1285905e6a2fSBarry Smith 12865615d1e5SSatish Balay #undef __FUNC__ 1287b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCopy_SeqDense" 1288cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 12894b0e389bSBarry Smith { 12904b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 12913a40ed3dSBarry Smith int ierr; 12923a40ed3dSBarry Smith 12933a40ed3dSBarry Smith PetscFunctionBegin; 12943a40ed3dSBarry Smith if (B->type != MATSEQDENSE) { 1295cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 12963a40ed3dSBarry Smith PetscFunctionReturn(0); 12973a40ed3dSBarry Smith } 1298*29bbc08cSBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1299549d3d68SSatish Balay ierr = PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));CHKERRQ(ierr); 13003a40ed3dSBarry Smith PetscFunctionReturn(0); 13014b0e389bSBarry Smith } 13024b0e389bSBarry Smith 1303289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1304a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1305905e6a2fSBarry Smith MatGetRow_SeqDense, 1306905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1307905e6a2fSBarry Smith MatMult_SeqDense, 1308905e6a2fSBarry Smith MatMultAdd_SeqDense, 13097c922b88SBarry Smith MatMultTranspose_SeqDense, 13107c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1311905e6a2fSBarry Smith MatSolve_SeqDense, 1312905e6a2fSBarry Smith MatSolveAdd_SeqDense, 13137c922b88SBarry Smith MatSolveTranspose_SeqDense, 13147c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1315905e6a2fSBarry Smith MatLUFactor_SeqDense, 1316905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1317ec8511deSBarry Smith MatRelax_SeqDense, 1318ec8511deSBarry Smith MatTranspose_SeqDense, 1319905e6a2fSBarry Smith MatGetInfo_SeqDense, 1320905e6a2fSBarry Smith MatEqual_SeqDense, 1321905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1322905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1323905e6a2fSBarry Smith MatNorm_SeqDense, 1324905e6a2fSBarry Smith 0, 1325905e6a2fSBarry Smith 0, 1326905e6a2fSBarry Smith 0, 1327905e6a2fSBarry Smith MatSetOption_SeqDense, 1328905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1329905e6a2fSBarry Smith MatZeroRows_SeqDense, 1330905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1331905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1332905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1333905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1334905e6a2fSBarry Smith MatGetSize_SeqDense, 1335905e6a2fSBarry Smith MatGetSize_SeqDense, 1336905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1337905e6a2fSBarry Smith 0, 1338905e6a2fSBarry Smith 0, 1339905e6a2fSBarry Smith MatGetArray_SeqDense, 1340905e6a2fSBarry Smith MatRestoreArray_SeqDense, 13415609ef8eSBarry Smith MatDuplicate_SeqDense, 1342a5ae1ecdSBarry Smith 0, 1343a5ae1ecdSBarry Smith 0, 1344a5ae1ecdSBarry Smith 0, 1345a5ae1ecdSBarry Smith 0, 1346a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1347a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1348a5ae1ecdSBarry Smith 0, 13494b0e389bSBarry Smith MatGetValues_SeqDense, 1350a5ae1ecdSBarry Smith MatCopy_SeqDense, 1351a5ae1ecdSBarry Smith 0, 1352a5ae1ecdSBarry Smith MatScale_SeqDense, 1353a5ae1ecdSBarry Smith 0, 1354a5ae1ecdSBarry Smith 0, 1355a5ae1ecdSBarry Smith 0, 1356a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1357a5ae1ecdSBarry Smith 0, 1358a5ae1ecdSBarry Smith 0, 1359a5ae1ecdSBarry Smith 0, 1360a5ae1ecdSBarry Smith 0, 1361a5ae1ecdSBarry Smith 0, 1362a5ae1ecdSBarry Smith 0, 1363a5ae1ecdSBarry Smith 0, 1364a5ae1ecdSBarry Smith 0, 1365a5ae1ecdSBarry Smith 0, 1366a5ae1ecdSBarry Smith 0, 1367e03a110bSBarry Smith MatDestroy_SeqDense, 1368e03a110bSBarry Smith MatView_SeqDense, 1369a5ae1ecdSBarry Smith MatGetMaps_Petsc}; 137090ace30eSBarry Smith 13715615d1e5SSatish Balay #undef __FUNC__ 1372b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateSeqDense" 13734b828684SBarry Smith /*@C 1374fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1375d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1376d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1377289bc588SBarry Smith 1378db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1379db81eaa0SLois Curfman McInnes 138020563c6bSBarry Smith Input Parameters: 1381db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 13820c775827SLois Curfman McInnes . m - number of rows 138318f449edSLois Curfman McInnes . n - number of columns 1384db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1385dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 138620563c6bSBarry Smith 138720563c6bSBarry Smith Output Parameter: 138844cd7ae7SLois Curfman McInnes . A - the matrix 138920563c6bSBarry Smith 1390b259b22eSLois Curfman McInnes Notes: 139118f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 139218f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1393b4fd4287SBarry Smith set data=PETSC_NULL. 139418f449edSLois Curfman McInnes 1395027ccd11SLois Curfman McInnes Level: intermediate 1396027ccd11SLois Curfman McInnes 1397dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1398d65003e9SLois Curfman McInnes 1399db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 140020563c6bSBarry Smith @*/ 140144cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1402289bc588SBarry Smith { 140344cd7ae7SLois Curfman McInnes Mat B; 140444cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 1405f1af5d2fSBarry Smith int ierr,size; 1406f1af5d2fSBarry Smith PetscTruth flg; 14073b2fbd54SBarry Smith 14083a40ed3dSBarry Smith PetscFunctionBegin; 1409d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1410*29bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 141155659b69SBarry Smith 141244cd7ae7SLois Curfman McInnes *A = 0; 14133f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,"Mat",comm,MatDestroy,MatView); 141444cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 141544cd7ae7SLois Curfman McInnes b = (Mat_SeqDense*)PetscMalloc(sizeof(Mat_SeqDense));CHKPTRQ(b); 1416549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1417549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 141844cd7ae7SLois Curfman McInnes B->factor = 0; 141990f02eecSBarry Smith B->mapping = 0; 1420f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 142144cd7ae7SLois Curfman McInnes B->data = (void*)b; 142218f449edSLois Curfman McInnes 142344cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 142444cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1425a5ae1ecdSBarry Smith 1426488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 1427488ecbafSBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1428a5ae1ecdSBarry Smith 142944cd7ae7SLois Curfman McInnes b->pivots = 0; 143044cd7ae7SLois Curfman McInnes b->roworiented = 1; 14317c922b88SBarry Smith if (!data) { 143244cd7ae7SLois Curfman McInnes b->v = (Scalar*)PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(b->v); 1433549d3d68SSatish Balay ierr = PetscMemzero(b->v,m*n*sizeof(Scalar));CHKERRQ(ierr); 143444cd7ae7SLois Curfman McInnes b->user_alloc = 0; 143544cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 1436273fdc76SBarry Smith } else { /* user-allocated storage */ 143744cd7ae7SLois Curfman McInnes b->v = data; 143844cd7ae7SLois Curfman McInnes b->user_alloc = 1; 14392dd2b441SLois Curfman McInnes } 144025cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 144144cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr);} 14424e220ebcSLois Curfman McInnes 144344cd7ae7SLois Curfman McInnes *A = B; 14443a40ed3dSBarry Smith PetscFunctionReturn(0); 1445289bc588SBarry Smith } 1446289bc588SBarry Smith 14473b2fbd54SBarry Smith 14483b2fbd54SBarry Smith 14493b2fbd54SBarry Smith 14503b2fbd54SBarry Smith 14513b2fbd54SBarry Smith 1452