1*ffac6cdbSBarry Smith /*$Id: dense.c,v 1.190 2000/10/24 20:25:29 bsmith 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; 14273d9f13SBarry 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; 27273d9f13SBarry Smith 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 33273d9f13SBarry Smith info->rows_global = (double)A->m; 34273d9f13SBarry Smith info->columns_global = (double)A->n; 35273d9f13SBarry Smith info->rows_local = (double)A->m; 36273d9f13SBarry Smith 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" 53273d9f13SBarry Smith int MatScale_SeqDense(Scalar *alpha,Mat A) 5480cd9d93SLois Curfman McInnes { 55273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 5680cd9d93SLois Curfman McInnes int one = 1,nz; 5780cd9d93SLois Curfman McInnes 583a40ed3dSBarry Smith PetscFunctionBegin; 59273d9f13SBarry Smith 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) { 77273d9f13SBarry Smith mat->pivots = (int*)PetscMalloc((A->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 78273d9f13SBarry Smith PLogObjectMemory(A,A->m*sizeof(int)); 79289bc588SBarry Smith } 80f1af5d2fSBarry Smith A->factor = FACTOR_LU; 81273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 82273d9f13SBarry Smith LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info); 8329bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 8429bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 85273d9f13SBarry Smith PLogFlops((2*A->n*A->n*A->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; 98273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr); 9902cad45dSBarry Smith l = (Mat_SeqDense*)newi->data; 1005609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 101273d9f13SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->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 */ 128273d9f13SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->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); 166273d9f13SBarry Smith PLogObjectMemory(A,-A->m*sizeof(int)); 167752f5567SLois Curfman McInnes mat->pivots = 0; 168752f5567SLois Curfman McInnes } 169273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 170273d9f13SBarry Smith LApotrf_("L",&A->n,mat->v,&A->m,&info); 17129bbc08cSBarry 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; 173273d9f13SBarry Smith PLogFlops((A->n*A->n*A->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; 186273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 1877a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1887a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 189273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 190c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 191273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 1927a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 193273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 194289bc588SBarry Smith } 19529bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 19629bbc08cSBarry 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); 199273d9f13SBarry Smith PLogFlops(2*A->n*A->n - A->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; 212273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2137a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2147a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 215273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 216752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 217da3a660dSBarry Smith if (mat->pivots) { 218273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 2197a97a34bSBarry Smith } else { 220273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 221da3a660dSBarry Smith } 22229bbc08cSBarry 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); 225273d9f13SBarry Smith PLogFlops(2*A->n*A->n - A->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); 241273d9f13SBarry Smith if (!A->m || !A->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 } 247273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 248752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 249da3a660dSBarry Smith if (mat->pivots) { 250273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 251a8c6a408SBarry Smith } else { 252273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 253da3a660dSBarry Smith } 25429bbc08cSBarry 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); 259273d9f13SBarry Smith PLogFlops(2*A->n*A->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; 273273d9f13SBarry Smith if (!A->m || !A->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 } 281273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 282752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 283da3a660dSBarry Smith if (mat->pivots) { 284273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 2853a40ed3dSBarry Smith } else { 286273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 287da3a660dSBarry Smith } 28829bbc08cSBarry 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); 297273d9f13SBarry Smith PLogFlops(2*A->n*A->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; 308273d9f13SBarry Smith int ierr,m = A->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; 371273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3727a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3737a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 374273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 3757a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3767a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 377273d9f13SBarry Smith PLogFlops(2*A->m*A->n - A->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; 390273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3917a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3927a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 393273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 3947a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3957a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 396273d9f13SBarry Smith PLogFlops(2*A->m*A->n - A->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; 409273d9f13SBarry Smith if (!A->m || !A->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); 413273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4147a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4157a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 416273d9f13SBarry Smith PLogFlops(2*A->m*A->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; 430273d9f13SBarry Smith if (!A->m || !A->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); 434273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4357a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4367a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 437273d9f13SBarry Smith PLogFlops(2*A->m*A->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; 451273d9f13SBarry Smith *ncols = A->n; 452289bc588SBarry Smith if (cols) { 453273d9f13SBarry Smith *cols = (int*)PetscMalloc((A->n+1)*sizeof(int));CHKPTRQ(*cols); 454273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 455289bc588SBarry Smith } 456289bc588SBarry Smith if (vals) { 457273d9f13SBarry Smith *vals = (Scalar*)PetscMalloc((A->n+1)*sizeof(Scalar));CHKPTRQ(*vals); 458289bc588SBarry Smith v = mat->v + row; 459273d9f13SBarry Smith for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->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) 48929bbc08cSBarry 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) 49429bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 49558804f6dSBarry Smith #endif 496273d9f13SBarry Smith mat->v[indexn[j]*A->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) 50329bbc08cSBarry 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) 50829bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 50958804f6dSBarry Smith #endif 510273d9f13SBarry Smith mat->v[indexn[j]*A->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) 51929bbc08cSBarry 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) 52429bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 52558804f6dSBarry Smith #endif 526273d9f13SBarry Smith mat->v[indexn[j]*A->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) 53329bbc08cSBarry 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) 53829bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 53958804f6dSBarry Smith #endif 540273d9f13SBarry Smith mat->v[indexn[j]*A->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++) { 560273d9f13SBarry Smith *vpt++ = mat->v[indexn[j]*A->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 570273d9f13SBarry Smith EXTERN_C_BEGIN 5715615d1e5SSatish Balay #undef __FUNC__ 572b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_SeqDense" 57319bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 57488e32edaSLois Curfman McInnes { 57588e32edaSLois Curfman McInnes Mat_SeqDense *a; 57688e32edaSLois Curfman McInnes Mat B; 57788e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 57888e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 5794a41a4d3SSatish Balay Scalar *vals,*svals,*v,*w; 58019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 58188e32edaSLois Curfman McInnes 5823a40ed3dSBarry Smith PetscFunctionBegin; 583d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 58429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 58590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 5860752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 58729bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 58888e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 58988e32edaSLois Curfman McInnes 59090ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 59190ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 59290ace30eSBarry Smith B = *A; 59390ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 5944a41a4d3SSatish Balay v = a->v; 5954a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 5964a41a4d3SSatish Balay from row major to column major */ 5974a41a4d3SSatish Balay w = (Scalar*)PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w); 59890ace30eSBarry Smith /* read in nonzero values */ 5994a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6004a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6014a41a4d3SSatish Balay for (j=0; j<N; j++) { 6024a41a4d3SSatish Balay for (i=0; i<M; i++) { 6034a41a4d3SSatish Balay *v++ =w[i*N+j]; 6044a41a4d3SSatish Balay } 6054a41a4d3SSatish Balay } 606606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6076d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6086d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60990ace30eSBarry Smith } else { 61088e32edaSLois Curfman McInnes /* read row lengths */ 61145d48df9SBarry Smith rowlengths = (int*)PetscMalloc((M+1)*sizeof(int));CHKPTRQ(rowlengths); 6120752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 61388e32edaSLois Curfman McInnes 61488e32edaSLois Curfman McInnes /* create our matrix */ 615b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 61688e32edaSLois Curfman McInnes B = *A; 61788e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 61888e32edaSLois Curfman McInnes v = a->v; 61988e32edaSLois Curfman McInnes 62088e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 62145d48df9SBarry Smith cols = scols = (int*)PetscMalloc((nz+1)*sizeof(int));CHKPTRQ(cols); 6220752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 62345d48df9SBarry Smith vals = svals = (Scalar*)PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(vals); 6240752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 62588e32edaSLois Curfman McInnes 62688e32edaSLois Curfman McInnes /* insert into matrix */ 62788e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 62888e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 62988e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 63088e32edaSLois Curfman McInnes } 631606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 632606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 633606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 63488e32edaSLois Curfman McInnes 6356d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6366d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63790ace30eSBarry Smith } 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 63988e32edaSLois Curfman McInnes } 640273d9f13SBarry Smith EXTERN_C_END 64188e32edaSLois Curfman McInnes 642e090d566SSatish Balay #include "petscsys.h" 643932b0c3eSLois Curfman McInnes 6445615d1e5SSatish Balay #undef __FUNC__ 645b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_ASCII" 646932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 647289bc588SBarry Smith { 648932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 649932b0c3eSLois Curfman McInnes int ierr,i,j,format; 650932b0c3eSLois Curfman McInnes char *outputname; 651932b0c3eSLois Curfman McInnes Scalar *v; 652932b0c3eSLois Curfman McInnes 6533a40ed3dSBarry Smith PetscFunctionBegin; 65477ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 655888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 656639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 6573a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 6586831982aSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 6596831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 660273d9f13SBarry Smith for (i=0; i<A->m; i++) { 66144cd7ae7SLois Curfman McInnes v = a->v + i; 6626831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 663273d9f13SBarry Smith for (j=0; j<A->n; j++) { 664aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 665329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 666329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 667329f5518SBarry Smith } else if (PetscRealPart(*v)) { 668329f5518SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));CHKERRQ(ierr); 6696831982aSBarry Smith } 67080cd9d93SLois Curfman McInnes #else 6716831982aSBarry Smith if (*v) { 6726831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr); 6736831982aSBarry Smith } 67480cd9d93SLois Curfman McInnes #endif 675273d9f13SBarry Smith v += A->m; 67680cd9d93SLois Curfman McInnes } 6776831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 67880cd9d93SLois Curfman McInnes } 6796831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 6803a40ed3dSBarry Smith } else { 6816831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 682aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 683*ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 68447989497SBarry Smith /* determine if matrix has all real values */ 68547989497SBarry Smith v = a->v; 686273d9f13SBarry Smith for (i=0; i<A->m*A->n; i++) { 687*ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 68847989497SBarry Smith } 68947989497SBarry Smith #endif 690*ffac6cdbSBarry Smith if (format == VIEWER_FORMAT_ASCII_MATLAB) { 691*ffac6cdbSBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 692*ffac6cdbSBarry Smith if (!outputname) outputname = "A"; 693*ffac6cdbSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 694*ffac6cdbSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",outputname,A->m,A->n);CHKERRQ(ierr); 695*ffac6cdbSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%s = [\n",outputname);CHKERRQ(ierr); 696*ffac6cdbSBarry Smith } 697*ffac6cdbSBarry Smith 698273d9f13SBarry Smith for (i=0; i<A->m; i++) { 699932b0c3eSLois Curfman McInnes v = a->v + i; 700273d9f13SBarry Smith for (j=0; j<A->n; j++) { 701aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 70247989497SBarry Smith if (allreal) { 703273d9f13SBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m; 70447989497SBarry Smith } else { 705273d9f13SBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m; 70647989497SBarry Smith } 707289bc588SBarry Smith #else 708273d9f13SBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m; 709289bc588SBarry Smith #endif 710289bc588SBarry Smith } 7116831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 712289bc588SBarry Smith } 713*ffac6cdbSBarry Smith if (format == VIEWER_FORMAT_ASCII_MATLAB) { 714*ffac6cdbSBarry Smith ierr = ViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 715*ffac6cdbSBarry Smith } 7166831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 717da3a660dSBarry Smith } 7186831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 7193a40ed3dSBarry Smith PetscFunctionReturn(0); 720289bc588SBarry Smith } 721289bc588SBarry Smith 7225615d1e5SSatish Balay #undef __FUNC__ 723b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Binary" 724932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 725932b0c3eSLois Curfman McInnes { 726932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 727273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 72890ace30eSBarry Smith int format; 72990ace30eSBarry Smith Scalar *v,*anonz,*vals; 730932b0c3eSLois Curfman McInnes 7313a40ed3dSBarry Smith PetscFunctionBegin; 73290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 73390ace30eSBarry Smith 73490ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 73502cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 73690ace30eSBarry Smith /* store the matrix as a dense matrix */ 73790ace30eSBarry Smith col_lens = (int*)PetscMalloc(4*sizeof(int));CHKPTRQ(col_lens); 73890ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 73990ace30eSBarry Smith col_lens[1] = m; 74090ace30eSBarry Smith col_lens[2] = n; 74190ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 7420752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 743606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 74490ace30eSBarry Smith 74590ace30eSBarry Smith /* write out matrix, by rows */ 74645d48df9SBarry Smith vals = (Scalar*)PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(vals); 74790ace30eSBarry Smith v = a->v; 74890ace30eSBarry Smith for (i=0; i<m; i++) { 74990ace30eSBarry Smith for (j=0; j<n; j++) { 75090ace30eSBarry Smith vals[i + j*m] = *v++; 75190ace30eSBarry Smith } 75290ace30eSBarry Smith } 7530752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 754606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 75590ace30eSBarry Smith } else { 7560452661fSBarry Smith col_lens = (int*)PetscMalloc((4+nz)*sizeof(int));CHKPTRQ(col_lens); 757932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 758932b0c3eSLois Curfman McInnes col_lens[1] = m; 759932b0c3eSLois Curfman McInnes col_lens[2] = n; 760932b0c3eSLois Curfman McInnes col_lens[3] = nz; 761932b0c3eSLois Curfman McInnes 762932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 763932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 7640752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 765932b0c3eSLois Curfman McInnes 766932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 767932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 768932b0c3eSLois Curfman McInnes ict = 0; 769932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 770932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 771932b0c3eSLois Curfman McInnes } 7720752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 773606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 774932b0c3eSLois Curfman McInnes 775932b0c3eSLois Curfman McInnes /* store nonzero values */ 77645d48df9SBarry Smith anonz = (Scalar*)PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(anonz); 777932b0c3eSLois Curfman McInnes ict = 0; 778932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 779932b0c3eSLois Curfman McInnes v = a->v + i; 780932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 781273d9f13SBarry Smith anonz[ict++] = *v; v += A->m; 782932b0c3eSLois Curfman McInnes } 783932b0c3eSLois Curfman McInnes } 7840752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 785606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 78690ace30eSBarry Smith } 7873a40ed3dSBarry Smith PetscFunctionReturn(0); 788932b0c3eSLois Curfman McInnes } 789932b0c3eSLois Curfman McInnes 7905615d1e5SSatish Balay #undef __FUNC__ 791b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Draw_Zoom" 792f1af5d2fSBarry Smith int MatView_SeqDense_Draw_Zoom(Draw draw,void *Aa) 793f1af5d2fSBarry Smith { 794f1af5d2fSBarry Smith Mat A = (Mat) Aa; 795f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 796273d9f13SBarry Smith int m = A->m,n = A->n,format,color,i,j,ierr; 797f1af5d2fSBarry Smith Scalar *v = a->v; 798f1af5d2fSBarry Smith Viewer viewer; 799f1af5d2fSBarry Smith Draw popup; 800329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 801f1af5d2fSBarry Smith 802f1af5d2fSBarry Smith PetscFunctionBegin; 803f1af5d2fSBarry Smith 804f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 805f1af5d2fSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 806f1af5d2fSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 807f1af5d2fSBarry Smith 808f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 809f1af5d2fSBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 810f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 811f1af5d2fSBarry Smith color = DRAW_BLUE; 812f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 813f1af5d2fSBarry Smith x_l = j; 814f1af5d2fSBarry Smith x_r = x_l + 1.0; 815f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 816f1af5d2fSBarry Smith y_l = m - i - 1.0; 817f1af5d2fSBarry Smith y_r = y_l + 1.0; 818f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 819329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 820f1af5d2fSBarry Smith color = DRAW_RED; 821329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 822f1af5d2fSBarry Smith color = DRAW_BLUE; 823f1af5d2fSBarry Smith } else { 824f1af5d2fSBarry Smith continue; 825f1af5d2fSBarry Smith } 826f1af5d2fSBarry Smith #else 827f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 828f1af5d2fSBarry Smith color = DRAW_RED; 829f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 830f1af5d2fSBarry Smith color = DRAW_BLUE; 831f1af5d2fSBarry Smith } else { 832f1af5d2fSBarry Smith continue; 833f1af5d2fSBarry Smith } 834f1af5d2fSBarry Smith #endif 835f1af5d2fSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 836f1af5d2fSBarry Smith } 837f1af5d2fSBarry Smith } 838f1af5d2fSBarry Smith } else { 839f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 840f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 841f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 842f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 843f1af5d2fSBarry Smith } 844f1af5d2fSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 845f1af5d2fSBarry Smith ierr = DrawGetPopup(draw,&popup);CHKERRQ(ierr); 8467c922b88SBarry Smith if (popup) {ierr = DrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 847f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 848f1af5d2fSBarry Smith x_l = j; 849f1af5d2fSBarry Smith x_r = x_l + 1.0; 850f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 851f1af5d2fSBarry Smith y_l = m - i - 1.0; 852f1af5d2fSBarry Smith y_r = y_l + 1.0; 853f1af5d2fSBarry Smith color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 854f1af5d2fSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 855f1af5d2fSBarry Smith } 856f1af5d2fSBarry Smith } 857f1af5d2fSBarry Smith } 858f1af5d2fSBarry Smith PetscFunctionReturn(0); 859f1af5d2fSBarry Smith } 860f1af5d2fSBarry Smith 861f1af5d2fSBarry Smith #undef __FUNC__ 862b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense_Draw" 863f1af5d2fSBarry Smith int MatView_SeqDense_Draw(Mat A,Viewer viewer) 864f1af5d2fSBarry Smith { 865f1af5d2fSBarry Smith Draw draw; 866f1af5d2fSBarry Smith PetscTruth isnull; 867329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 868f1af5d2fSBarry Smith int ierr; 869f1af5d2fSBarry Smith 870f1af5d2fSBarry Smith PetscFunctionBegin; 871f1af5d2fSBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 872f1af5d2fSBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); 873f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 874f1af5d2fSBarry Smith 875f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 876273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 877f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 878f1af5d2fSBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 879f1af5d2fSBarry Smith ierr = DrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 880f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 881f1af5d2fSBarry Smith PetscFunctionReturn(0); 882f1af5d2fSBarry Smith } 883f1af5d2fSBarry Smith 884f1af5d2fSBarry Smith #undef __FUNC__ 885b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqDense" 886c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer) 887932b0c3eSLois Curfman McInnes { 888932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 889bcd2baecSBarry Smith int ierr; 890f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 891932b0c3eSLois Curfman McInnes 8923a40ed3dSBarry Smith PetscFunctionBegin; 8936831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 8946831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 8956831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 896f1af5d2fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 8970f5bd95cSBarry Smith 8980f5bd95cSBarry Smith if (issocket) { 899273d9f13SBarry Smith ierr = ViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9000f5bd95cSBarry Smith } else if (isascii) { 9013a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9020f5bd95cSBarry Smith } else if (isbinary) { 9033a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 904f1af5d2fSBarry Smith } else if (isdraw) { 905f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9065cd90555SBarry Smith } else { 90729bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 908932b0c3eSLois Curfman McInnes } 9093a40ed3dSBarry Smith PetscFunctionReturn(0); 910932b0c3eSLois Curfman McInnes } 911289bc588SBarry Smith 9125615d1e5SSatish Balay #undef __FUNC__ 913b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_SeqDense" 914c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 915289bc588SBarry Smith { 916ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 91790f02eecSBarry Smith int ierr; 91890f02eecSBarry Smith 9193a40ed3dSBarry Smith PetscFunctionBegin; 920aa482453SBarry Smith #if defined(PETSC_USE_LOG) 921273d9f13SBarry Smith PLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 922a5a9c739SBarry Smith #endif 923606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 924606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 925606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 9263a40ed3dSBarry Smith PetscFunctionReturn(0); 927289bc588SBarry Smith } 928289bc588SBarry Smith 9295615d1e5SSatish Balay #undef __FUNC__ 930b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_SeqDense" 9318f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 932289bc588SBarry Smith { 933c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 934549d3d68SSatish Balay int k,j,m,n,ierr; 935d3e5ee88SLois Curfman McInnes Scalar *v,tmp; 93648b35521SBarry Smith 9373a40ed3dSBarry Smith PetscFunctionBegin; 938273d9f13SBarry Smith v = mat->v; m = A->m; n = A->n; 9397c922b88SBarry Smith if (!matout) { /* in place transpose */ 940d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 941d64ed03dSBarry Smith Scalar *w = (Scalar*)PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 942d64ed03dSBarry Smith for (j=0; j<m; j++) { 943d64ed03dSBarry Smith for (k=0; k<n; k++) { 944d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 945d64ed03dSBarry Smith } 946d64ed03dSBarry Smith } 947549d3d68SSatish Balay ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr); 948606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 949d64ed03dSBarry Smith } else { 950d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 951289bc588SBarry Smith for (k=0; k<j; k++) { 952d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 953d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 954d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 955289bc588SBarry Smith } 956289bc588SBarry Smith } 957d64ed03dSBarry Smith } 9583a40ed3dSBarry Smith } else { /* out-of-place transpose */ 959d3e5ee88SLois Curfman McInnes Mat tmat; 960ec8511deSBarry Smith Mat_SeqDense *tmatd; 961d3e5ee88SLois Curfman McInnes Scalar *v2; 962273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 963ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 9640de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 965d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 9660de55854SLois Curfman McInnes for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m]; 967d3e5ee88SLois Curfman McInnes } 9686d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9696d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 970d3e5ee88SLois Curfman McInnes *matout = tmat; 97148b35521SBarry Smith } 9723a40ed3dSBarry Smith PetscFunctionReturn(0); 973289bc588SBarry Smith } 974289bc588SBarry Smith 9755615d1e5SSatish Balay #undef __FUNC__ 976b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatEqual_SeqDense" 9778f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 978289bc588SBarry Smith { 979c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 980c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 981273d9f13SBarry Smith int ierr,i; 982289bc588SBarry Smith Scalar *v1 = mat1->v,*v2 = mat2->v; 983273d9f13SBarry Smith PetscTruth flag; 9849ea5d5aeSSatish Balay 9853a40ed3dSBarry Smith PetscFunctionBegin; 986273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr); 987273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 988273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 989273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 990273d9f13SBarry Smith for (i=0; i<A1->m*A1->n; i++) { 9913a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 992289bc588SBarry Smith v1++; v2++; 993289bc588SBarry Smith } 99477c4ece6SBarry Smith *flg = PETSC_TRUE; 9953a40ed3dSBarry Smith PetscFunctionReturn(0); 996289bc588SBarry Smith } 997289bc588SBarry Smith 9985615d1e5SSatish Balay #undef __FUNC__ 999b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetDiagonal_SeqDense" 10008f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1001289bc588SBarry Smith { 1002c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10037a97a34bSBarry Smith int ierr,i,n,len; 100444cd7ae7SLois Curfman McInnes Scalar *x,zero = 0.0; 100544cd7ae7SLois Curfman McInnes 10063a40ed3dSBarry Smith PetscFunctionBegin; 10077a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10087a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1009600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1010273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1011273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 101244cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 1013273d9f13SBarry Smith x[i] = mat->v[i*A->m + i]; 1014289bc588SBarry Smith } 10157a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10163a40ed3dSBarry Smith PetscFunctionReturn(0); 1017289bc588SBarry Smith } 1018289bc588SBarry Smith 10195615d1e5SSatish Balay #undef __FUNC__ 1020b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDiagonalScale_SeqDense" 10218f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1022289bc588SBarry Smith { 1023c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1024da3a660dSBarry Smith Scalar *l,*r,x,*v; 1025273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 102655659b69SBarry Smith 10273a40ed3dSBarry Smith PetscFunctionBegin; 102828988994SBarry Smith if (ll) { 10297a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1030600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1031273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1032da3a660dSBarry Smith for (i=0; i<m; i++) { 1033da3a660dSBarry Smith x = l[i]; 1034da3a660dSBarry Smith v = mat->v + i; 1035da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1036da3a660dSBarry Smith } 10377a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 103844cd7ae7SLois Curfman McInnes PLogFlops(n*m); 1039da3a660dSBarry Smith } 104028988994SBarry Smith if (rr) { 10417a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1042600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1043273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1044da3a660dSBarry Smith for (i=0; i<n; i++) { 1045da3a660dSBarry Smith x = r[i]; 1046da3a660dSBarry Smith v = mat->v + i*m; 1047da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1048da3a660dSBarry Smith } 10497a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 105044cd7ae7SLois Curfman McInnes PLogFlops(n*m); 1051da3a660dSBarry Smith } 10523a40ed3dSBarry Smith PetscFunctionReturn(0); 1053289bc588SBarry Smith } 1054289bc588SBarry Smith 10555615d1e5SSatish Balay #undef __FUNC__ 1056b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatNorm_SeqDense" 1057329f5518SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm) 1058289bc588SBarry Smith { 1059c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1060289bc588SBarry Smith Scalar *v = mat->v; 1061329f5518SBarry Smith PetscReal sum = 0.0; 1062289bc588SBarry Smith int i,j; 106355659b69SBarry Smith 10643a40ed3dSBarry Smith PetscFunctionBegin; 1065289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1066273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1067aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1068329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1069289bc588SBarry Smith #else 1070289bc588SBarry Smith sum += (*v)*(*v); v++; 1071289bc588SBarry Smith #endif 1072289bc588SBarry Smith } 1073289bc588SBarry Smith *norm = sqrt(sum); 1074273d9f13SBarry Smith PLogFlops(2*A->n*A->m); 10753a40ed3dSBarry Smith } else if (type == NORM_1) { 1076289bc588SBarry Smith *norm = 0.0; 1077273d9f13SBarry Smith for (j=0; j<A->n; j++) { 1078289bc588SBarry Smith sum = 0.0; 1079273d9f13SBarry Smith for (i=0; i<A->m; i++) { 108033a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1081289bc588SBarry Smith } 1082289bc588SBarry Smith if (sum > *norm) *norm = sum; 1083289bc588SBarry Smith } 1084273d9f13SBarry Smith PLogFlops(A->n*A->m); 10853a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1086289bc588SBarry Smith *norm = 0.0; 1087273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1088289bc588SBarry Smith v = mat->v + j; 1089289bc588SBarry Smith sum = 0.0; 1090273d9f13SBarry Smith for (i=0; i<A->n; i++) { 1091273d9f13SBarry Smith sum += PetscAbsScalar(*v); v += A->m; 1092289bc588SBarry Smith } 1093289bc588SBarry Smith if (sum > *norm) *norm = sum; 1094289bc588SBarry Smith } 1095273d9f13SBarry Smith PLogFlops(A->n*A->m); 10963a40ed3dSBarry Smith } else { 109729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1098289bc588SBarry Smith } 10993a40ed3dSBarry Smith PetscFunctionReturn(0); 1100289bc588SBarry Smith } 1101289bc588SBarry Smith 11025615d1e5SSatish Balay #undef __FUNC__ 1103b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_SeqDense" 11048f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1105289bc588SBarry Smith { 1106c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 110767e560aaSBarry Smith 11083a40ed3dSBarry Smith PetscFunctionBegin; 1109273d9f13SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = PETSC_TRUE; 1110273d9f13SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = PETSC_FALSE; 11116d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 1112219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 11136d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 1114219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 11156d4a8577SBarry Smith op == MAT_SYMMETRIC || 11166d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 11176d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 11186d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 11194787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 11206d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 112190f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 1122b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 1123c38d4ed2SBarry Smith op == MAT_USE_HASH_TABLE) { 1124981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1125c38d4ed2SBarry Smith } else { 112629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 11273a40ed3dSBarry Smith } 11283a40ed3dSBarry Smith PetscFunctionReturn(0); 1129289bc588SBarry Smith } 1130289bc588SBarry Smith 11315615d1e5SSatish Balay #undef __FUNC__ 1132b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroEntries_SeqDense" 11338f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 11346f0a148fSBarry Smith { 1135ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1136549d3d68SSatish Balay int ierr; 11373a40ed3dSBarry Smith 11383a40ed3dSBarry Smith PetscFunctionBegin; 1139273d9f13SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 11403a40ed3dSBarry Smith PetscFunctionReturn(0); 11416f0a148fSBarry Smith } 11426f0a148fSBarry Smith 11435615d1e5SSatish Balay #undef __FUNC__ 1144b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_SeqDense" 11458f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 11464e220ebcSLois Curfman McInnes { 11473a40ed3dSBarry Smith PetscFunctionBegin; 11484e220ebcSLois Curfman McInnes *bs = 1; 11493a40ed3dSBarry Smith PetscFunctionReturn(0); 11504e220ebcSLois Curfman McInnes } 11514e220ebcSLois Curfman McInnes 11525615d1e5SSatish Balay #undef __FUNC__ 1153b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqDense" 11548f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 11556f0a148fSBarry Smith { 1156ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1157273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 11586f0a148fSBarry Smith Scalar *slot; 115955659b69SBarry Smith 11603a40ed3dSBarry Smith PetscFunctionBegin; 1161e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 116278b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 11636f0a148fSBarry Smith for (i=0; i<N; i++) { 11646f0a148fSBarry Smith slot = l->v + rows[i]; 11656f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 11666f0a148fSBarry Smith } 11676f0a148fSBarry Smith if (diag) { 11686f0a148fSBarry Smith for (i=0; i<N; i++) { 11696f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1170260925b8SBarry Smith *slot = *diag; 11716f0a148fSBarry Smith } 11726f0a148fSBarry Smith } 1173260925b8SBarry Smith ISRestoreIndices(is,&rows); 11743a40ed3dSBarry Smith PetscFunctionReturn(0); 11756f0a148fSBarry Smith } 1176557bce09SLois Curfman McInnes 11775615d1e5SSatish Balay #undef __FUNC__ 1178b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_SeqDense" 11798f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1180d3e5ee88SLois Curfman McInnes { 11813a40ed3dSBarry Smith PetscFunctionBegin; 11824c49b128SBarry Smith if (m) *m = 0; 1183273d9f13SBarry Smith if (n) *n = A->m; 11843a40ed3dSBarry Smith PetscFunctionReturn(0); 1185d3e5ee88SLois Curfman McInnes } 1186d3e5ee88SLois Curfman McInnes 11875615d1e5SSatish Balay #undef __FUNC__ 1188b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetArray_SeqDense" 11898f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 119064e87e97SBarry Smith { 1191c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11923a40ed3dSBarry Smith 11933a40ed3dSBarry Smith PetscFunctionBegin; 119464e87e97SBarry Smith *array = mat->v; 11953a40ed3dSBarry Smith PetscFunctionReturn(0); 119664e87e97SBarry Smith } 11970754003eSLois Curfman McInnes 11985615d1e5SSatish Balay #undef __FUNC__ 1199b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreArray_SeqDense" 12008f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1201ff14e315SSatish Balay { 12023a40ed3dSBarry Smith PetscFunctionBegin; 120309b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12043a40ed3dSBarry Smith PetscFunctionReturn(0); 1205ff14e315SSatish Balay } 12060754003eSLois Curfman McInnes 12075615d1e5SSatish Balay #undef __FUNC__ 1208b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrix_SeqDense" 12097b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12100754003eSLois Curfman McInnes { 1211c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1212273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 1213182d2002SSatish Balay Scalar *av,*bv,*v = mat->v; 12140754003eSLois Curfman McInnes Mat newmat; 12150754003eSLois Curfman McInnes 12163a40ed3dSBarry Smith PetscFunctionBegin; 121778b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 121878b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1219e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1220e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 12210754003eSLois Curfman McInnes 1222182d2002SSatish Balay /* Check submatrixcall */ 1223182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1224182d2002SSatish Balay int n_cols,n_rows; 1225182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 122629bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1227182d2002SSatish Balay newmat = *B; 1228182d2002SSatish Balay } else { 12290754003eSLois Curfman McInnes /* Create and fill new matrix */ 1230b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1231182d2002SSatish Balay } 1232182d2002SSatish Balay 1233182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1234182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1235182d2002SSatish Balay 1236182d2002SSatish Balay for (i=0; i<ncols; i++) { 1237182d2002SSatish Balay av = v + m*icol[i]; 1238182d2002SSatish Balay for (j=0; j<nrows; j++) { 1239182d2002SSatish Balay *bv++ = av[irow[j]]; 12400754003eSLois Curfman McInnes } 12410754003eSLois Curfman McInnes } 1242182d2002SSatish Balay 1243182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 12446d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12456d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12460754003eSLois Curfman McInnes 12470754003eSLois Curfman McInnes /* Free work space */ 124878b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 124978b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1250182d2002SSatish Balay *B = newmat; 12513a40ed3dSBarry Smith PetscFunctionReturn(0); 12520754003eSLois Curfman McInnes } 12530754003eSLois Curfman McInnes 12545615d1e5SSatish Balay #undef __FUNC__ 1255b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSubMatrices_SeqDense" 12567b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1257905e6a2fSBarry Smith { 1258905e6a2fSBarry Smith int ierr,i; 1259905e6a2fSBarry Smith 12603a40ed3dSBarry Smith PetscFunctionBegin; 1261905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1262905e6a2fSBarry Smith *B = (Mat*)PetscMalloc((n+1)*sizeof(Mat));CHKPTRQ(*B); 1263905e6a2fSBarry Smith } 1264905e6a2fSBarry Smith 1265905e6a2fSBarry Smith for (i=0; i<n; i++) { 12666a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1267905e6a2fSBarry Smith } 12683a40ed3dSBarry Smith PetscFunctionReturn(0); 1269905e6a2fSBarry Smith } 1270905e6a2fSBarry Smith 12715615d1e5SSatish Balay #undef __FUNC__ 1272b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCopy_SeqDense" 1273cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 12744b0e389bSBarry Smith { 12754b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 12763a40ed3dSBarry Smith int ierr; 1277273d9f13SBarry Smith PetscTruth flag; 12783a40ed3dSBarry Smith 12793a40ed3dSBarry Smith PetscFunctionBegin; 1280273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr); 1281273d9f13SBarry Smith if (!flag) { 1282cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 12833a40ed3dSBarry Smith PetscFunctionReturn(0); 12843a40ed3dSBarry Smith } 1285273d9f13SBarry Smith if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1286273d9f13SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 1287273d9f13SBarry Smith PetscFunctionReturn(0); 1288273d9f13SBarry Smith } 1289273d9f13SBarry Smith 1290273d9f13SBarry Smith #undef __FUNC__ 1291273d9f13SBarry Smith #define __FUNC__ /*<a name="MatSetUpPreallocation_SeqDense"></a>*/"MatSetUpPreallocation_SeqDense" 1292273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1293273d9f13SBarry Smith { 1294273d9f13SBarry Smith int ierr; 1295273d9f13SBarry Smith 1296273d9f13SBarry Smith PetscFunctionBegin; 1297273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 12983a40ed3dSBarry Smith PetscFunctionReturn(0); 12994b0e389bSBarry Smith } 13004b0e389bSBarry Smith 1301289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1302a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1303905e6a2fSBarry Smith MatGetRow_SeqDense, 1304905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1305905e6a2fSBarry Smith MatMult_SeqDense, 1306905e6a2fSBarry Smith MatMultAdd_SeqDense, 13077c922b88SBarry Smith MatMultTranspose_SeqDense, 13087c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1309905e6a2fSBarry Smith MatSolve_SeqDense, 1310905e6a2fSBarry Smith MatSolveAdd_SeqDense, 13117c922b88SBarry Smith MatSolveTranspose_SeqDense, 13127c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1313905e6a2fSBarry Smith MatLUFactor_SeqDense, 1314905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1315ec8511deSBarry Smith MatRelax_SeqDense, 1316ec8511deSBarry Smith MatTranspose_SeqDense, 1317905e6a2fSBarry Smith MatGetInfo_SeqDense, 1318905e6a2fSBarry Smith MatEqual_SeqDense, 1319905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1320905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1321905e6a2fSBarry Smith MatNorm_SeqDense, 1322905e6a2fSBarry Smith 0, 1323905e6a2fSBarry Smith 0, 1324905e6a2fSBarry Smith 0, 1325905e6a2fSBarry Smith MatSetOption_SeqDense, 1326905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1327905e6a2fSBarry Smith MatZeroRows_SeqDense, 1328905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1329905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1330905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1331905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1332273d9f13SBarry Smith MatSetUpPreallocation_SeqDense, 1333273d9f13SBarry Smith 0, 1334905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1335905e6a2fSBarry Smith 0, 1336905e6a2fSBarry Smith 0, 1337905e6a2fSBarry Smith MatGetArray_SeqDense, 1338905e6a2fSBarry Smith MatRestoreArray_SeqDense, 13395609ef8eSBarry Smith MatDuplicate_SeqDense, 1340a5ae1ecdSBarry Smith 0, 1341a5ae1ecdSBarry Smith 0, 1342a5ae1ecdSBarry Smith 0, 1343a5ae1ecdSBarry Smith 0, 1344a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1345a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1346a5ae1ecdSBarry Smith 0, 13474b0e389bSBarry Smith MatGetValues_SeqDense, 1348a5ae1ecdSBarry Smith MatCopy_SeqDense, 1349a5ae1ecdSBarry Smith 0, 1350a5ae1ecdSBarry Smith MatScale_SeqDense, 1351a5ae1ecdSBarry Smith 0, 1352a5ae1ecdSBarry Smith 0, 1353a5ae1ecdSBarry Smith 0, 1354a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1355a5ae1ecdSBarry Smith 0, 1356a5ae1ecdSBarry Smith 0, 1357a5ae1ecdSBarry Smith 0, 1358a5ae1ecdSBarry Smith 0, 1359a5ae1ecdSBarry Smith 0, 1360a5ae1ecdSBarry Smith 0, 1361a5ae1ecdSBarry Smith 0, 1362a5ae1ecdSBarry Smith 0, 1363a5ae1ecdSBarry Smith 0, 1364a5ae1ecdSBarry Smith 0, 1365e03a110bSBarry Smith MatDestroy_SeqDense, 1366e03a110bSBarry Smith MatView_SeqDense, 1367a5ae1ecdSBarry Smith MatGetMaps_Petsc}; 136890ace30eSBarry Smith 13695615d1e5SSatish Balay #undef __FUNC__ 1370b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateSeqDense" 13714b828684SBarry Smith /*@C 1372fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1373d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1374d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1375289bc588SBarry Smith 1376db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1377db81eaa0SLois Curfman McInnes 137820563c6bSBarry Smith Input Parameters: 1379db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 13800c775827SLois Curfman McInnes . m - number of rows 138118f449edSLois Curfman McInnes . n - number of columns 1382db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1383dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 138420563c6bSBarry Smith 138520563c6bSBarry Smith Output Parameter: 138644cd7ae7SLois Curfman McInnes . A - the matrix 138720563c6bSBarry Smith 1388b259b22eSLois Curfman McInnes Notes: 138918f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 139018f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1391b4fd4287SBarry Smith set data=PETSC_NULL. 139218f449edSLois Curfman McInnes 1393027ccd11SLois Curfman McInnes Level: intermediate 1394027ccd11SLois Curfman McInnes 1395dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1396d65003e9SLois Curfman McInnes 1397db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 139820563c6bSBarry Smith @*/ 139944cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1400289bc588SBarry Smith { 1401273d9f13SBarry Smith int ierr; 14023b2fbd54SBarry Smith 14033a40ed3dSBarry Smith PetscFunctionBegin; 1404273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1405273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1406273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1407273d9f13SBarry Smith PetscFunctionReturn(0); 1408273d9f13SBarry Smith } 1409273d9f13SBarry Smith 1410273d9f13SBarry Smith #undef __FUNC__ 1411273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqDensePreallocation" 1412273d9f13SBarry Smith /*@C 1413273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1414273d9f13SBarry Smith 1415273d9f13SBarry Smith Collective on MPI_Comm 1416273d9f13SBarry Smith 1417273d9f13SBarry Smith Input Parameters: 1418273d9f13SBarry Smith + A - the matrix 1419273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1420273d9f13SBarry Smith 1421273d9f13SBarry Smith Notes: 1422273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1423273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1424273d9f13SBarry Smith set data=PETSC_NULL. 1425273d9f13SBarry Smith 1426273d9f13SBarry Smith Level: intermediate 1427273d9f13SBarry Smith 1428273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1429273d9f13SBarry Smith 1430273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1431273d9f13SBarry Smith @*/ 1432273d9f13SBarry Smith int MatSeqDenseSetPreallocation(Mat B,Scalar *data) 1433273d9f13SBarry Smith { 1434273d9f13SBarry Smith Mat_SeqDense *b; 1435273d9f13SBarry Smith int ierr; 1436273d9f13SBarry Smith PetscTruth flg2; 1437273d9f13SBarry Smith 1438273d9f13SBarry Smith PetscFunctionBegin; 1439273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr); 1440273d9f13SBarry Smith if (!flg2) PetscFunctionReturn(0); 1441273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1442273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1443273d9f13SBarry Smith if (!data) { 1444273d9f13SBarry Smith b->v = (Scalar*)PetscMalloc((B->m*B->n+1)*sizeof(Scalar));CHKPTRQ(b->v); 1445273d9f13SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(Scalar));CHKERRQ(ierr); 1446273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 1447273d9f13SBarry Smith PLogObjectMemory(B,B->n*B->m*sizeof(Scalar)); 1448273d9f13SBarry Smith } else { /* user-allocated storage */ 1449273d9f13SBarry Smith b->v = data; 1450273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1451273d9f13SBarry Smith } 1452273d9f13SBarry Smith PetscFunctionReturn(0); 1453273d9f13SBarry Smith } 1454273d9f13SBarry Smith 1455273d9f13SBarry Smith EXTERN_C_BEGIN 1456273d9f13SBarry Smith #undef __FUNC__ 1457273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_SeqDense" 1458273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1459273d9f13SBarry Smith { 1460273d9f13SBarry Smith Mat_SeqDense *b; 1461273d9f13SBarry Smith int ierr,size; 1462273d9f13SBarry Smith 1463273d9f13SBarry Smith PetscFunctionBegin; 1464273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 146529bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 146655659b69SBarry Smith 1467273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1468273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1469273d9f13SBarry Smith 147044cd7ae7SLois Curfman McInnes b = (Mat_SeqDense*)PetscMalloc(sizeof(Mat_SeqDense));CHKPTRQ(b); 1471549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1472549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 147344cd7ae7SLois Curfman McInnes B->factor = 0; 147490f02eecSBarry Smith B->mapping = 0; 1475f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 147644cd7ae7SLois Curfman McInnes B->data = (void*)b; 147718f449edSLois Curfman McInnes 1478273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 1479273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1480a5ae1ecdSBarry Smith 148144cd7ae7SLois Curfman McInnes b->pivots = 0; 1482273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1483273d9f13SBarry Smith b->v = 0; 14844e220ebcSLois Curfman McInnes 14853a40ed3dSBarry Smith PetscFunctionReturn(0); 1486289bc588SBarry Smith } 1487273d9f13SBarry Smith EXTERN_C_END 14883b2fbd54SBarry Smith 14893b2fbd54SBarry Smith 14903b2fbd54SBarry Smith 14913b2fbd54SBarry Smith 1492