1*ae7cfcebSSatish Balay /*$Id: dense.c,v 1.197 2001/03/23 23:21:48 balay Exp balay $*/ 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 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "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); 18b0a32e0cSBarry Smith PetscLogFlops(2*N-1); 193a40ed3dSBarry Smith PetscFunctionReturn(0); 201987afe7SBarry Smith } 211987afe7SBarry Smith 224a2ae208SSatish Balay #undef __FUNCT__ 234a2ae208SSatish Balay #define __FUNCT__ "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 514a2ae208SSatish Balay #undef __FUNCT__ 524a2ae208SSatish Balay #define __FUNCT__ "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); 61b0a32e0cSBarry Smith PetscLogFlops(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.*/ 684a2ae208SSatish Balay #undef __FUNCT__ 694a2ae208SSatish Balay #define __FUNCT__ "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; 73b0a32e0cSBarry Smith int info,ierr; 7455659b69SBarry Smith 753a40ed3dSBarry Smith PetscFunctionBegin; 76289bc588SBarry Smith if (!mat->pivots) { 77b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr); 78b0a32e0cSBarry Smith PetscLogObjectMemory(A,A->m*sizeof(int)); 79289bc588SBarry Smith } 80f1af5d2fSBarry Smith A->factor = FACTOR_LU; 81273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 82*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRF) 83*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable."); 84*ae7cfcebSSatish Balay #else 85273d9f13SBarry Smith LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info); 86*ae7cfcebSSatish Balay #endif 8729bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 8829bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 89b0a32e0cSBarry Smith PetscLogFlops((2*A->n*A->n*A->n)/3); 903a40ed3dSBarry Smith PetscFunctionReturn(0); 91289bc588SBarry Smith } 926ee01492SSatish Balay 934a2ae208SSatish Balay #undef __FUNCT__ 944a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 955609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9602cad45dSBarry Smith { 9702cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 9802cad45dSBarry Smith int ierr; 9902cad45dSBarry Smith Mat newi; 10002cad45dSBarry Smith 1013a40ed3dSBarry Smith PetscFunctionBegin; 102273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr); 10302cad45dSBarry Smith l = (Mat_SeqDense*)newi->data; 1045609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 105273d9f13SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 10602cad45dSBarry Smith } 10702cad45dSBarry Smith newi->assembled = PETSC_TRUE; 10802cad45dSBarry Smith *newmat = newi; 1093a40ed3dSBarry Smith PetscFunctionReturn(0); 11002cad45dSBarry Smith } 11102cad45dSBarry Smith 1124a2ae208SSatish Balay #undef __FUNCT__ 1134a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 114e03a110bSBarry Smith int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatLUInfo *info,Mat *fact) 115289bc588SBarry Smith { 1163a40ed3dSBarry Smith int ierr; 1173a40ed3dSBarry Smith 1183a40ed3dSBarry Smith PetscFunctionBegin; 1195609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1203a40ed3dSBarry Smith PetscFunctionReturn(0); 121289bc588SBarry Smith } 1226ee01492SSatish Balay 1234a2ae208SSatish Balay #undef __FUNCT__ 1244a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1258f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 126289bc588SBarry Smith { 12702cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 1283a40ed3dSBarry Smith int ierr; 1293a40ed3dSBarry Smith 1303a40ed3dSBarry Smith PetscFunctionBegin; 13102cad45dSBarry Smith /* copy the numerical values */ 132273d9f13SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 13302cad45dSBarry Smith (*fact)->factor = 0; 134e03a110bSBarry Smith ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr); 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 136289bc588SBarry Smith } 1376ee01492SSatish Balay 1384a2ae208SSatish Balay #undef __FUNCT__ 1394a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 140329f5518SBarry Smith int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,PetscReal f,Mat *fact) 141289bc588SBarry Smith { 1423a40ed3dSBarry Smith int ierr; 1433a40ed3dSBarry Smith 1443a40ed3dSBarry Smith PetscFunctionBegin; 1453a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1463a40ed3dSBarry Smith PetscFunctionReturn(0); 147289bc588SBarry Smith } 1486ee01492SSatish Balay 1494a2ae208SSatish Balay #undef __FUNCT__ 1504a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 1518f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 152289bc588SBarry Smith { 1533a40ed3dSBarry Smith int ierr; 1543a40ed3dSBarry Smith 1553a40ed3dSBarry Smith PetscFunctionBegin; 1563a40ed3dSBarry Smith ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr); 1573a40ed3dSBarry Smith PetscFunctionReturn(0); 158289bc588SBarry Smith } 1596ee01492SSatish Balay 1604a2ae208SSatish Balay #undef __FUNCT__ 1614a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 162329f5518SBarry Smith int MatCholeskyFactor_SeqDense(Mat A,IS perm,PetscReal f) 163289bc588SBarry Smith { 164c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 165606d414cSSatish Balay int info,ierr; 16655659b69SBarry Smith 1673a40ed3dSBarry Smith PetscFunctionBegin; 168752f5567SLois Curfman McInnes if (mat->pivots) { 169606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 170b0a32e0cSBarry Smith PetscLogObjectMemory(A,-A->m*sizeof(int)); 171752f5567SLois Curfman McInnes mat->pivots = 0; 172752f5567SLois Curfman McInnes } 173273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 174*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRF) 175*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable."); 176*ae7cfcebSSatish Balay #else 177273d9f13SBarry Smith LApotrf_("L",&A->n,mat->v,&A->m,&info); 178*ae7cfcebSSatish Balay #endif 17929bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 180c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 181b0a32e0cSBarry Smith PetscLogFlops((A->n*A->n*A->n)/3); 1823a40ed3dSBarry Smith PetscFunctionReturn(0); 183289bc588SBarry Smith } 184289bc588SBarry Smith 1854a2ae208SSatish Balay #undef __FUNCT__ 1864a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 1878f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 188289bc588SBarry Smith { 189c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 190a57ad546SLois Curfman McInnes int one = 1,info,ierr; 1916abc6512SBarry Smith Scalar *x,*y; 19267e560aaSBarry Smith 1933a40ed3dSBarry Smith PetscFunctionBegin; 194273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 1957a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1967a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 197273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 198c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 199*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 200*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 201*ae7cfcebSSatish Balay #else 202273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 203*ae7cfcebSSatish Balay #endif 2047a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 205*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 206*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 207*ae7cfcebSSatish Balay #else 208273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 209*ae7cfcebSSatish Balay #endif 210289bc588SBarry Smith } 21129bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 21229bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 2137a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2147a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 215b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2163a40ed3dSBarry Smith PetscFunctionReturn(0); 217289bc588SBarry Smith } 2186ee01492SSatish Balay 2194a2ae208SSatish Balay #undef __FUNCT__ 2204a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 2217c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 222da3a660dSBarry Smith { 223c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2247a97a34bSBarry Smith int ierr,one = 1,info; 2256abc6512SBarry Smith Scalar *x,*y; 22667e560aaSBarry Smith 2273a40ed3dSBarry Smith PetscFunctionBegin; 228273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2297a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2307a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 231273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 232752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 233da3a660dSBarry Smith if (mat->pivots) { 234*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 235*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 236*ae7cfcebSSatish Balay #else 237273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 238*ae7cfcebSSatish Balay #endif 2397a97a34bSBarry Smith } else { 240*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 241*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 242*ae7cfcebSSatish Balay #else 243273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 244*ae7cfcebSSatish Balay #endif 245da3a660dSBarry Smith } 24629bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 2477a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2487a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 249b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2503a40ed3dSBarry Smith PetscFunctionReturn(0); 251da3a660dSBarry Smith } 2526ee01492SSatish Balay 2534a2ae208SSatish Balay #undef __FUNCT__ 2544a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 2558f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 256da3a660dSBarry Smith { 257c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2587a97a34bSBarry Smith int ierr,one = 1,info; 2596abc6512SBarry Smith Scalar *x,*y,sone = 1.0; 260da3a660dSBarry Smith Vec tmp = 0; 26167e560aaSBarry Smith 2623a40ed3dSBarry Smith PetscFunctionBegin; 2637a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2647a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 265273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 266da3a660dSBarry Smith if (yy == zz) { 26778b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 268b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 26978b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 270da3a660dSBarry Smith } 271273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 272752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 273da3a660dSBarry Smith if (mat->pivots) { 274*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 275*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 276*ae7cfcebSSatish Balay #else 277273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 278*ae7cfcebSSatish Balay #endif 279a8c6a408SBarry Smith } else { 280*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 281*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 282*ae7cfcebSSatish Balay #else 283273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 284*ae7cfcebSSatish Balay #endif 285da3a660dSBarry Smith } 28629bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 287a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 288a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 2897a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2907a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 291b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 2923a40ed3dSBarry Smith PetscFunctionReturn(0); 293da3a660dSBarry Smith } 29467e560aaSBarry Smith 2954a2ae208SSatish Balay #undef __FUNCT__ 2964a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 2977c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 298da3a660dSBarry Smith { 299c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3006abc6512SBarry Smith int one = 1,info,ierr; 3016abc6512SBarry Smith Scalar *x,*y,sone = 1.0; 302da3a660dSBarry Smith Vec tmp; 30367e560aaSBarry Smith 3043a40ed3dSBarry Smith PetscFunctionBegin; 305273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3067a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3077a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 308da3a660dSBarry Smith if (yy == zz) { 30978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 310b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 31178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 312da3a660dSBarry Smith } 313273d9f13SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(Scalar));CHKERRQ(ierr); 314752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 315da3a660dSBarry Smith if (mat->pivots) { 316*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 317*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 318*ae7cfcebSSatish Balay #else 319273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 320*ae7cfcebSSatish Balay #endif 3213a40ed3dSBarry Smith } else { 322*ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 323*ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 324*ae7cfcebSSatish Balay #else 325273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 326*ae7cfcebSSatish Balay #endif 327da3a660dSBarry Smith } 32829bbc08cSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 32990f02eecSBarry Smith if (tmp) { 33090f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 33190f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3323a40ed3dSBarry Smith } else { 33390f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 33490f02eecSBarry Smith } 3357a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3367a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 337b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 339da3a660dSBarry Smith } 340289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3414a2ae208SSatish Balay #undef __FUNCT__ 3424a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 343329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 344329f5518SBarry Smith PetscReal shift,int its,Vec xx) 345289bc588SBarry Smith { 346c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 347289bc588SBarry Smith Scalar *x,*b,*v = mat->v,zero = 0.0,xt; 348273d9f13SBarry Smith int ierr,m = A->m,i; 349aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 350bc1b551cSSatish Balay int o = 1; 351bc1b551cSSatish Balay #endif 352289bc588SBarry Smith 3533a40ed3dSBarry Smith PetscFunctionBegin; 354289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3553a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 356bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 357289bc588SBarry Smith } 3587a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3597a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 360289bc588SBarry Smith while (its--) { 361289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 362289bc588SBarry Smith for (i=0; i<m; i++) { 363aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 364f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 365f1747703SBarry Smith not happy about returning a double complex */ 366f1747703SBarry Smith int _i; 367f1747703SBarry Smith Scalar sum = b[i]; 368f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3693f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 370f1747703SBarry Smith } 371f1747703SBarry Smith xt = sum; 372f1747703SBarry Smith #else 373289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 374f1747703SBarry Smith #endif 37555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 376289bc588SBarry Smith } 377289bc588SBarry Smith } 378289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 379289bc588SBarry Smith for (i=m-1; i>=0; i--) { 380aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 381f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 382f1747703SBarry Smith not happy about returning a double complex */ 383f1747703SBarry Smith int _i; 384f1747703SBarry Smith Scalar sum = b[i]; 385f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3863f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 387f1747703SBarry Smith } 388f1747703SBarry Smith xt = sum; 389f1747703SBarry Smith #else 390289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 391f1747703SBarry Smith #endif 39255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 393289bc588SBarry Smith } 394289bc588SBarry Smith } 395289bc588SBarry Smith } 396600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 3977a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3983a40ed3dSBarry Smith PetscFunctionReturn(0); 399289bc588SBarry Smith } 400289bc588SBarry Smith 401289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4024a2ae208SSatish Balay #undef __FUNCT__ 4034a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 4047c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 405289bc588SBarry Smith { 406c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 407289bc588SBarry Smith Scalar *v = mat->v,*x,*y; 4087a97a34bSBarry Smith int ierr,_One=1;Scalar _DOne=1.0,_DZero=0.0; 4093a40ed3dSBarry Smith 4103a40ed3dSBarry Smith PetscFunctionBegin; 411273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4127a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4137a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 414273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 4157a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4167a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 417b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->n); 4183a40ed3dSBarry Smith PetscFunctionReturn(0); 419289bc588SBarry Smith } 4206ee01492SSatish Balay 4214a2ae208SSatish Balay #undef __FUNCT__ 4224a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 42344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 424289bc588SBarry Smith { 425c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 426329f5518SBarry Smith Scalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 427329f5518SBarry Smith int ierr,_One=1; 4283a40ed3dSBarry Smith 4293a40ed3dSBarry Smith PetscFunctionBegin; 430273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4317a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4327a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 433273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 4347a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4357a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 436b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->m); 4373a40ed3dSBarry Smith PetscFunctionReturn(0); 438289bc588SBarry Smith } 4396ee01492SSatish Balay 4404a2ae208SSatish Balay #undef __FUNCT__ 4414a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 44244cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 443289bc588SBarry Smith { 444c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 445329f5518SBarry Smith Scalar *v = mat->v,*x,*y,_DOne=1.0; 446329f5518SBarry Smith int ierr,_One=1; 4473a40ed3dSBarry Smith 4483a40ed3dSBarry Smith PetscFunctionBegin; 449273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 450600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4517a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4527a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 453273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4547a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4557a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 456b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4573a40ed3dSBarry Smith PetscFunctionReturn(0); 458289bc588SBarry Smith } 4596ee01492SSatish Balay 4604a2ae208SSatish Balay #undef __FUNCT__ 4614a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 4627c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 463289bc588SBarry Smith { 464c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 465600d6d0bSBarry Smith Scalar *v = mat->v,*x,*y; 4667a97a34bSBarry Smith int ierr,_One=1; 4676abc6512SBarry Smith Scalar _DOne=1.0; 4683a40ed3dSBarry Smith 4693a40ed3dSBarry Smith PetscFunctionBegin; 470273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 471600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4727a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4737a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 474273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4757a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4767a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 477b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4783a40ed3dSBarry Smith PetscFunctionReturn(0); 479289bc588SBarry Smith } 480289bc588SBarry Smith 481289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4824a2ae208SSatish Balay #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 4848f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 485289bc588SBarry Smith { 486c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 487289bc588SBarry Smith Scalar *v; 488b0a32e0cSBarry Smith int i,ierr; 48967e560aaSBarry Smith 4903a40ed3dSBarry Smith PetscFunctionBegin; 491273d9f13SBarry Smith *ncols = A->n; 492289bc588SBarry Smith if (cols) { 493b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 494273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 495289bc588SBarry Smith } 496289bc588SBarry Smith if (vals) { 497b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(Scalar),vals);CHKERRQ(ierr); 498289bc588SBarry Smith v = mat->v + row; 499273d9f13SBarry Smith for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;} 500289bc588SBarry Smith } 5013a40ed3dSBarry Smith PetscFunctionReturn(0); 502289bc588SBarry Smith } 5036ee01492SSatish Balay 5044a2ae208SSatish Balay #undef __FUNCT__ 5054a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 5068f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 507289bc588SBarry Smith { 508606d414cSSatish Balay int ierr; 509606d414cSSatish Balay PetscFunctionBegin; 510606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 511606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5123a40ed3dSBarry Smith PetscFunctionReturn(0); 513289bc588SBarry Smith } 514289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5154a2ae208SSatish Balay #undef __FUNCT__ 5164a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 5178f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 518e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 519289bc588SBarry Smith { 520c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 521289bc588SBarry Smith int i,j; 522d6dfbf8fSBarry Smith 5233a40ed3dSBarry Smith PetscFunctionBegin; 524289bc588SBarry Smith if (!mat->roworiented) { 525dbb450caSBarry Smith if (addv == INSERT_VALUES) { 526289bc588SBarry Smith for (j=0; j<n; j++) { 5275ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 528aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 52929bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 53058804f6dSBarry Smith #endif 531289bc588SBarry Smith for (i=0; i<m; i++) { 5325ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 533aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 53429bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 53558804f6dSBarry Smith #endif 536273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] = *v++; 537289bc588SBarry Smith } 538289bc588SBarry Smith } 5393a40ed3dSBarry Smith } else { 540289bc588SBarry Smith for (j=0; j<n; j++) { 5415ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 542aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 54329bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 54458804f6dSBarry Smith #endif 545289bc588SBarry Smith for (i=0; i<m; i++) { 5465ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 547aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 54829bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 54958804f6dSBarry Smith #endif 550273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] += *v++; 551289bc588SBarry Smith } 552289bc588SBarry Smith } 553289bc588SBarry Smith } 5543a40ed3dSBarry Smith } else { 555dbb450caSBarry Smith if (addv == INSERT_VALUES) { 556e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 5575ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 558aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 55929bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56058804f6dSBarry Smith #endif 561e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 5625ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 563aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 56429bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56558804f6dSBarry Smith #endif 566273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] = *v++; 567e8d4e0b9SBarry Smith } 568e8d4e0b9SBarry Smith } 5693a40ed3dSBarry Smith } else { 570289bc588SBarry Smith for (i=0; i<m; i++) { 5715ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 572aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57329bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 57458804f6dSBarry Smith #endif 575289bc588SBarry Smith for (j=0; j<n; j++) { 5765ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 577aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57829bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 57958804f6dSBarry Smith #endif 580273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] += *v++; 581289bc588SBarry Smith } 582289bc588SBarry Smith } 583289bc588SBarry Smith } 584e8d4e0b9SBarry Smith } 5853a40ed3dSBarry Smith PetscFunctionReturn(0); 586289bc588SBarry Smith } 587e8d4e0b9SBarry Smith 5884a2ae208SSatish Balay #undef __FUNCT__ 5894a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 5908f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 591ae80bb75SLois Curfman McInnes { 592ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 593ae80bb75SLois Curfman McInnes int i,j; 594ae80bb75SLois Curfman McInnes Scalar *vpt = v; 595ae80bb75SLois Curfman McInnes 5963a40ed3dSBarry Smith PetscFunctionBegin; 597ae80bb75SLois Curfman McInnes /* row-oriented output */ 598ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 599ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 600273d9f13SBarry Smith *vpt++ = mat->v[indexn[j]*A->m + indexm[i]]; 601ae80bb75SLois Curfman McInnes } 602ae80bb75SLois Curfman McInnes } 6033a40ed3dSBarry Smith PetscFunctionReturn(0); 604ae80bb75SLois Curfman McInnes } 605ae80bb75SLois Curfman McInnes 606289bc588SBarry Smith /* -----------------------------------------------------------------*/ 607289bc588SBarry Smith 608e090d566SSatish Balay #include "petscsys.h" 60988e32edaSLois Curfman McInnes 610273d9f13SBarry Smith EXTERN_C_BEGIN 6114a2ae208SSatish Balay #undef __FUNCT__ 6124a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 613b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A) 61488e32edaSLois Curfman McInnes { 61588e32edaSLois Curfman McInnes Mat_SeqDense *a; 61688e32edaSLois Curfman McInnes Mat B; 61788e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 61888e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 6194a41a4d3SSatish Balay Scalar *vals,*svals,*v,*w; 62019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 62188e32edaSLois Curfman McInnes 6223a40ed3dSBarry Smith PetscFunctionBegin; 623d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 62429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 625b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6260752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 62729bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 62888e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 62988e32edaSLois Curfman McInnes 63090ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 63190ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 63290ace30eSBarry Smith B = *A; 63390ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6344a41a4d3SSatish Balay v = a->v; 6354a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6364a41a4d3SSatish Balay from row major to column major */ 637b0a32e0cSBarry Smith ierr = PetscMalloc((M*N+1)*sizeof(Scalar),&w);CHKERRQ(ierr); 63890ace30eSBarry Smith /* read in nonzero values */ 6394a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6404a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6414a41a4d3SSatish Balay for (j=0; j<N; j++) { 6424a41a4d3SSatish Balay for (i=0; i<M; i++) { 6434a41a4d3SSatish Balay *v++ =w[i*N+j]; 6444a41a4d3SSatish Balay } 6454a41a4d3SSatish Balay } 646606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6476d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6486d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64990ace30eSBarry Smith } else { 65088e32edaSLois Curfman McInnes /* read row lengths */ 651b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 6520752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 65388e32edaSLois Curfman McInnes 65488e32edaSLois Curfman McInnes /* create our matrix */ 655b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 65688e32edaSLois Curfman McInnes B = *A; 65788e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 65888e32edaSLois Curfman McInnes v = a->v; 65988e32edaSLois Curfman McInnes 66088e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 661b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 662b0a32e0cSBarry Smith cols = scols; 6630752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 664b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(Scalar),&svals);CHKERRQ(ierr); 665b0a32e0cSBarry Smith vals = svals; 6660752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 66788e32edaSLois Curfman McInnes 66888e32edaSLois Curfman McInnes /* insert into matrix */ 66988e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 67088e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 67188e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 67288e32edaSLois Curfman McInnes } 673606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 674606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 675606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 67688e32edaSLois Curfman McInnes 6776d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6786d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67990ace30eSBarry Smith } 6803a40ed3dSBarry Smith PetscFunctionReturn(0); 68188e32edaSLois Curfman McInnes } 682273d9f13SBarry Smith EXTERN_C_END 68388e32edaSLois Curfman McInnes 684e090d566SSatish Balay #include "petscsys.h" 685932b0c3eSLois Curfman McInnes 6864a2ae208SSatish Balay #undef __FUNCT__ 6874a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 688b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 689289bc588SBarry Smith { 690932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 691fb9695e5SSatish Balay int ierr,i,j; 692fb9695e5SSatish Balay char *name; 693932b0c3eSLois Curfman McInnes Scalar *v; 694f3ef73ceSBarry Smith PetscViewerFormat format; 695932b0c3eSLois Curfman McInnes 6963a40ed3dSBarry Smith PetscFunctionBegin; 697435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 698b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 699fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 7003a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 701fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 702b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 703273d9f13SBarry Smith for (i=0; i<A->m; i++) { 70444cd7ae7SLois Curfman McInnes v = a->v + i; 705b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 706273d9f13SBarry Smith for (j=0; j<A->n; j++) { 707aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 708329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 709b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 710329f5518SBarry Smith } else if (PetscRealPart(*v)) { 711b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,PetscRealPart(*v));CHKERRQ(ierr); 7126831982aSBarry Smith } 71380cd9d93SLois Curfman McInnes #else 7146831982aSBarry Smith if (*v) { 715b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr); 7166831982aSBarry Smith } 71780cd9d93SLois Curfman McInnes #endif 718273d9f13SBarry Smith v += A->m; 71980cd9d93SLois Curfman McInnes } 720b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 72180cd9d93SLois Curfman McInnes } 722b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7233a40ed3dSBarry Smith } else { 724b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 725aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 726ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 72747989497SBarry Smith /* determine if matrix has all real values */ 72847989497SBarry Smith v = a->v; 729273d9f13SBarry Smith for (i=0; i<A->m*A->n; i++) { 730ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 73147989497SBarry Smith } 73247989497SBarry Smith #endif 733fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 734fb9695e5SSatish Balay ierr = PetscObjectGetName((PetscObject)viewer,&name);CHKERRQ(ierr); 735b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 736fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 737fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 738ffac6cdbSBarry Smith } 739ffac6cdbSBarry Smith 740273d9f13SBarry Smith for (i=0; i<A->m; i++) { 741932b0c3eSLois Curfman McInnes v = a->v + i; 742273d9f13SBarry Smith for (j=0; j<A->n; j++) { 743aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 74447989497SBarry Smith if (allreal) { 745b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m; 74647989497SBarry Smith } else { 747b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m; 74847989497SBarry Smith } 749289bc588SBarry Smith #else 750b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m; 751289bc588SBarry Smith #endif 752289bc588SBarry Smith } 753b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 754289bc588SBarry Smith } 755fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 756b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 757ffac6cdbSBarry Smith } 758b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 759da3a660dSBarry Smith } 760b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7613a40ed3dSBarry Smith PetscFunctionReturn(0); 762289bc588SBarry Smith } 763289bc588SBarry Smith 7644a2ae208SSatish Balay #undef __FUNCT__ 7654a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 766b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 767932b0c3eSLois Curfman McInnes { 768932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 769273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 77090ace30eSBarry Smith Scalar *v,*anonz,*vals; 771f3ef73ceSBarry Smith PetscViewerFormat format; 772932b0c3eSLois Curfman McInnes 7733a40ed3dSBarry Smith PetscFunctionBegin; 774b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 77590ace30eSBarry Smith 776b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 777fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 77890ace30eSBarry Smith /* store the matrix as a dense matrix */ 779b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 78090ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 78190ace30eSBarry Smith col_lens[1] = m; 78290ace30eSBarry Smith col_lens[2] = n; 78390ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 7840752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 785606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 78690ace30eSBarry Smith 78790ace30eSBarry Smith /* write out matrix, by rows */ 788b0a32e0cSBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(Scalar),&vals);CHKERRQ(ierr); 78990ace30eSBarry Smith v = a->v; 79090ace30eSBarry Smith for (i=0; i<m; i++) { 79190ace30eSBarry Smith for (j=0; j<n; j++) { 79290ace30eSBarry Smith vals[i + j*m] = *v++; 79390ace30eSBarry Smith } 79490ace30eSBarry Smith } 7950752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 796606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 79790ace30eSBarry Smith } else { 798b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 799932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 800932b0c3eSLois Curfman McInnes col_lens[1] = m; 801932b0c3eSLois Curfman McInnes col_lens[2] = n; 802932b0c3eSLois Curfman McInnes col_lens[3] = nz; 803932b0c3eSLois Curfman McInnes 804932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 805932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8060752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 807932b0c3eSLois Curfman McInnes 808932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 809932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 810932b0c3eSLois Curfman McInnes ict = 0; 811932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 812932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 813932b0c3eSLois Curfman McInnes } 8140752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 815606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 816932b0c3eSLois Curfman McInnes 817932b0c3eSLois Curfman McInnes /* store nonzero values */ 818b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(Scalar),&anonz);CHKERRQ(ierr); 819932b0c3eSLois Curfman McInnes ict = 0; 820932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 821932b0c3eSLois Curfman McInnes v = a->v + i; 822932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 823273d9f13SBarry Smith anonz[ict++] = *v; v += A->m; 824932b0c3eSLois Curfman McInnes } 825932b0c3eSLois Curfman McInnes } 8260752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 827606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 82890ace30eSBarry Smith } 8293a40ed3dSBarry Smith PetscFunctionReturn(0); 830932b0c3eSLois Curfman McInnes } 831932b0c3eSLois Curfman McInnes 8324a2ae208SSatish Balay #undef __FUNCT__ 8334a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 834b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 835f1af5d2fSBarry Smith { 836f1af5d2fSBarry Smith Mat A = (Mat) Aa; 837f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 838fb9695e5SSatish Balay int m = A->m,n = A->n,color,i,j,ierr; 839f1af5d2fSBarry Smith Scalar *v = a->v; 840b0a32e0cSBarry Smith PetscViewer viewer; 841b0a32e0cSBarry Smith PetscDraw popup; 842329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 843f3ef73ceSBarry Smith PetscViewerFormat format; 844f1af5d2fSBarry Smith 845f1af5d2fSBarry Smith PetscFunctionBegin; 846f1af5d2fSBarry Smith 847f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 848b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 849b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 850f1af5d2fSBarry Smith 851f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 852fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 853f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 854b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 855f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 856f1af5d2fSBarry Smith x_l = j; 857f1af5d2fSBarry Smith x_r = x_l + 1.0; 858f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 859f1af5d2fSBarry Smith y_l = m - i - 1.0; 860f1af5d2fSBarry Smith y_r = y_l + 1.0; 861f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 862329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 863b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 864329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 865b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 866f1af5d2fSBarry Smith } else { 867f1af5d2fSBarry Smith continue; 868f1af5d2fSBarry Smith } 869f1af5d2fSBarry Smith #else 870f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 871b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 872f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 873b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 874f1af5d2fSBarry Smith } else { 875f1af5d2fSBarry Smith continue; 876f1af5d2fSBarry Smith } 877f1af5d2fSBarry Smith #endif 878b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 879f1af5d2fSBarry Smith } 880f1af5d2fSBarry Smith } 881f1af5d2fSBarry Smith } else { 882f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 883f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 884f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 885f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 886f1af5d2fSBarry Smith } 887b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 888b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 889b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 890f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 891f1af5d2fSBarry Smith x_l = j; 892f1af5d2fSBarry Smith x_r = x_l + 1.0; 893f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 894f1af5d2fSBarry Smith y_l = m - i - 1.0; 895f1af5d2fSBarry Smith y_r = y_l + 1.0; 896b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 897b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 898f1af5d2fSBarry Smith } 899f1af5d2fSBarry Smith } 900f1af5d2fSBarry Smith } 901f1af5d2fSBarry Smith PetscFunctionReturn(0); 902f1af5d2fSBarry Smith } 903f1af5d2fSBarry Smith 9044a2ae208SSatish Balay #undef __FUNCT__ 9054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 906b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 907f1af5d2fSBarry Smith { 908b0a32e0cSBarry Smith PetscDraw draw; 909f1af5d2fSBarry Smith PetscTruth isnull; 910329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 911f1af5d2fSBarry Smith int ierr; 912f1af5d2fSBarry Smith 913f1af5d2fSBarry Smith PetscFunctionBegin; 914b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 915b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 916f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 917f1af5d2fSBarry Smith 918f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 919273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 920f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 921b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 922b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 923f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 924f1af5d2fSBarry Smith PetscFunctionReturn(0); 925f1af5d2fSBarry Smith } 926f1af5d2fSBarry Smith 9274a2ae208SSatish Balay #undef __FUNCT__ 9284a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 929b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer) 930932b0c3eSLois Curfman McInnes { 931932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 932bcd2baecSBarry Smith int ierr; 933f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 934932b0c3eSLois Curfman McInnes 9353a40ed3dSBarry Smith PetscFunctionBegin; 936b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 937b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 938fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 939fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 9400f5bd95cSBarry Smith 9410f5bd95cSBarry Smith if (issocket) { 942b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9430f5bd95cSBarry Smith } else if (isascii) { 9443a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9450f5bd95cSBarry Smith } else if (isbinary) { 9463a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 947f1af5d2fSBarry Smith } else if (isdraw) { 948f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9495cd90555SBarry Smith } else { 95029bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 951932b0c3eSLois Curfman McInnes } 9523a40ed3dSBarry Smith PetscFunctionReturn(0); 953932b0c3eSLois Curfman McInnes } 954289bc588SBarry Smith 9554a2ae208SSatish Balay #undef __FUNCT__ 9564a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 957c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 958289bc588SBarry Smith { 959ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 96090f02eecSBarry Smith int ierr; 96190f02eecSBarry Smith 9623a40ed3dSBarry Smith PetscFunctionBegin; 963aa482453SBarry Smith #if defined(PETSC_USE_LOG) 964b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 965a5a9c739SBarry Smith #endif 966606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 967606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 968606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 9693a40ed3dSBarry Smith PetscFunctionReturn(0); 970289bc588SBarry Smith } 971289bc588SBarry Smith 9724a2ae208SSatish Balay #undef __FUNCT__ 9734a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 9748f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 975289bc588SBarry Smith { 976c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 977549d3d68SSatish Balay int k,j,m,n,ierr; 978d3e5ee88SLois Curfman McInnes Scalar *v,tmp; 97948b35521SBarry Smith 9803a40ed3dSBarry Smith PetscFunctionBegin; 981273d9f13SBarry Smith v = mat->v; m = A->m; n = A->n; 9827c922b88SBarry Smith if (!matout) { /* in place transpose */ 983d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 984b0a32e0cSBarry Smith Scalar *w; 985b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*(n+1)*sizeof(Scalar),&w);CHKERRQ(ierr); 986d64ed03dSBarry Smith for (j=0; j<m; j++) { 987d64ed03dSBarry Smith for (k=0; k<n; k++) { 988d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 989d64ed03dSBarry Smith } 990d64ed03dSBarry Smith } 991549d3d68SSatish Balay ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr); 992606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 993d64ed03dSBarry Smith } else { 994d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 995289bc588SBarry Smith for (k=0; k<j; k++) { 996d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 997d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 998d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 999289bc588SBarry Smith } 1000289bc588SBarry Smith } 1001d64ed03dSBarry Smith } 10023a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1003d3e5ee88SLois Curfman McInnes Mat tmat; 1004ec8511deSBarry Smith Mat_SeqDense *tmatd; 1005d3e5ee88SLois Curfman McInnes Scalar *v2; 1006273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1007ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10080de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1009d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10100de55854SLois Curfman McInnes for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m]; 1011d3e5ee88SLois Curfman McInnes } 10126d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10136d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1014d3e5ee88SLois Curfman McInnes *matout = tmat; 101548b35521SBarry Smith } 10163a40ed3dSBarry Smith PetscFunctionReturn(0); 1017289bc588SBarry Smith } 1018289bc588SBarry Smith 10194a2ae208SSatish Balay #undef __FUNCT__ 10204a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 10218f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1022289bc588SBarry Smith { 1023c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1024c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1025273d9f13SBarry Smith int ierr,i; 1026289bc588SBarry Smith Scalar *v1 = mat1->v,*v2 = mat2->v; 1027273d9f13SBarry Smith PetscTruth flag; 10289ea5d5aeSSatish Balay 10293a40ed3dSBarry Smith PetscFunctionBegin; 1030273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr); 1031273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 1032273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1033273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1034273d9f13SBarry Smith for (i=0; i<A1->m*A1->n; i++) { 10353a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1036289bc588SBarry Smith v1++; v2++; 1037289bc588SBarry Smith } 103877c4ece6SBarry Smith *flg = PETSC_TRUE; 10393a40ed3dSBarry Smith PetscFunctionReturn(0); 1040289bc588SBarry Smith } 1041289bc588SBarry Smith 10424a2ae208SSatish Balay #undef __FUNCT__ 10434a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 10448f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1045289bc588SBarry Smith { 1046c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10477a97a34bSBarry Smith int ierr,i,n,len; 104844cd7ae7SLois Curfman McInnes Scalar *x,zero = 0.0; 104944cd7ae7SLois Curfman McInnes 10503a40ed3dSBarry Smith PetscFunctionBegin; 10517a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10527a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1053600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1054273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1055273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 105644cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 1057273d9f13SBarry Smith x[i] = mat->v[i*A->m + i]; 1058289bc588SBarry Smith } 10597a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10603a40ed3dSBarry Smith PetscFunctionReturn(0); 1061289bc588SBarry Smith } 1062289bc588SBarry Smith 10634a2ae208SSatish Balay #undef __FUNCT__ 10644a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 10658f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1066289bc588SBarry Smith { 1067c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1068da3a660dSBarry Smith Scalar *l,*r,x,*v; 1069273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 107055659b69SBarry Smith 10713a40ed3dSBarry Smith PetscFunctionBegin; 107228988994SBarry Smith if (ll) { 10737a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1074600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1075273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1076da3a660dSBarry Smith for (i=0; i<m; i++) { 1077da3a660dSBarry Smith x = l[i]; 1078da3a660dSBarry Smith v = mat->v + i; 1079da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1080da3a660dSBarry Smith } 10817a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1082b0a32e0cSBarry Smith PetscLogFlops(n*m); 1083da3a660dSBarry Smith } 108428988994SBarry Smith if (rr) { 10857a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1086600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1087273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1088da3a660dSBarry Smith for (i=0; i<n; i++) { 1089da3a660dSBarry Smith x = r[i]; 1090da3a660dSBarry Smith v = mat->v + i*m; 1091da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1092da3a660dSBarry Smith } 10937a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1094b0a32e0cSBarry Smith PetscLogFlops(n*m); 1095da3a660dSBarry Smith } 10963a40ed3dSBarry Smith PetscFunctionReturn(0); 1097289bc588SBarry Smith } 1098289bc588SBarry Smith 10994a2ae208SSatish Balay #undef __FUNCT__ 11004a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1101329f5518SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *norm) 1102289bc588SBarry Smith { 1103c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1104289bc588SBarry Smith Scalar *v = mat->v; 1105329f5518SBarry Smith PetscReal sum = 0.0; 1106289bc588SBarry Smith int i,j; 110755659b69SBarry Smith 11083a40ed3dSBarry Smith PetscFunctionBegin; 1109289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1110273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1111aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1112329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1113289bc588SBarry Smith #else 1114289bc588SBarry Smith sum += (*v)*(*v); v++; 1115289bc588SBarry Smith #endif 1116289bc588SBarry Smith } 1117289bc588SBarry Smith *norm = sqrt(sum); 1118b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11193a40ed3dSBarry Smith } else if (type == NORM_1) { 1120289bc588SBarry Smith *norm = 0.0; 1121273d9f13SBarry Smith for (j=0; j<A->n; j++) { 1122289bc588SBarry Smith sum = 0.0; 1123273d9f13SBarry Smith for (i=0; i<A->m; i++) { 112433a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1125289bc588SBarry Smith } 1126289bc588SBarry Smith if (sum > *norm) *norm = sum; 1127289bc588SBarry Smith } 1128b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11293a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1130289bc588SBarry Smith *norm = 0.0; 1131273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1132289bc588SBarry Smith v = mat->v + j; 1133289bc588SBarry Smith sum = 0.0; 1134273d9f13SBarry Smith for (i=0; i<A->n; i++) { 1135273d9f13SBarry Smith sum += PetscAbsScalar(*v); v += A->m; 1136289bc588SBarry Smith } 1137289bc588SBarry Smith if (sum > *norm) *norm = sum; 1138289bc588SBarry Smith } 1139b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11403a40ed3dSBarry Smith } else { 114129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1142289bc588SBarry Smith } 11433a40ed3dSBarry Smith PetscFunctionReturn(0); 1144289bc588SBarry Smith } 1145289bc588SBarry Smith 11464a2ae208SSatish Balay #undef __FUNCT__ 11474a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 11488f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1149289bc588SBarry Smith { 1150c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 115167e560aaSBarry Smith 11523a40ed3dSBarry Smith PetscFunctionBegin; 1153273d9f13SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = PETSC_TRUE; 1154273d9f13SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = PETSC_FALSE; 11556d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 1156219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 11576d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 1158219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 11596d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 11606d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 11614787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 11626d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 116390f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 1164b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 1165c38d4ed2SBarry Smith op == MAT_USE_HASH_TABLE) { 1166b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1167c38d4ed2SBarry Smith } else { 116829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 11693a40ed3dSBarry Smith } 11703a40ed3dSBarry Smith PetscFunctionReturn(0); 1171289bc588SBarry Smith } 1172289bc588SBarry Smith 11734a2ae208SSatish Balay #undef __FUNCT__ 11744a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 11758f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 11766f0a148fSBarry Smith { 1177ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1178549d3d68SSatish Balay int ierr; 11793a40ed3dSBarry Smith 11803a40ed3dSBarry Smith PetscFunctionBegin; 1181273d9f13SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 11823a40ed3dSBarry Smith PetscFunctionReturn(0); 11836f0a148fSBarry Smith } 11846f0a148fSBarry Smith 11854a2ae208SSatish Balay #undef __FUNCT__ 11864a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 11878f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 11884e220ebcSLois Curfman McInnes { 11893a40ed3dSBarry Smith PetscFunctionBegin; 11904e220ebcSLois Curfman McInnes *bs = 1; 11913a40ed3dSBarry Smith PetscFunctionReturn(0); 11924e220ebcSLois Curfman McInnes } 11934e220ebcSLois Curfman McInnes 11944a2ae208SSatish Balay #undef __FUNCT__ 11954a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 11968f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 11976f0a148fSBarry Smith { 1198ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1199273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 12006f0a148fSBarry Smith Scalar *slot; 120155659b69SBarry Smith 12023a40ed3dSBarry Smith PetscFunctionBegin; 1203e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 120478b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12056f0a148fSBarry Smith for (i=0; i<N; i++) { 12066f0a148fSBarry Smith slot = l->v + rows[i]; 12076f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 12086f0a148fSBarry Smith } 12096f0a148fSBarry Smith if (diag) { 12106f0a148fSBarry Smith for (i=0; i<N; i++) { 12116f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1212260925b8SBarry Smith *slot = *diag; 12136f0a148fSBarry Smith } 12146f0a148fSBarry Smith } 1215260925b8SBarry Smith ISRestoreIndices(is,&rows); 12163a40ed3dSBarry Smith PetscFunctionReturn(0); 12176f0a148fSBarry Smith } 1218557bce09SLois Curfman McInnes 12194a2ae208SSatish Balay #undef __FUNCT__ 12204a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqDense" 12218f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1222d3e5ee88SLois Curfman McInnes { 12233a40ed3dSBarry Smith PetscFunctionBegin; 12244c49b128SBarry Smith if (m) *m = 0; 1225273d9f13SBarry Smith if (n) *n = A->m; 12263a40ed3dSBarry Smith PetscFunctionReturn(0); 1227d3e5ee88SLois Curfman McInnes } 1228d3e5ee88SLois Curfman McInnes 12294a2ae208SSatish Balay #undef __FUNCT__ 12304a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 12318f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 123264e87e97SBarry Smith { 1233c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12343a40ed3dSBarry Smith 12353a40ed3dSBarry Smith PetscFunctionBegin; 123664e87e97SBarry Smith *array = mat->v; 12373a40ed3dSBarry Smith PetscFunctionReturn(0); 123864e87e97SBarry Smith } 12390754003eSLois Curfman McInnes 12404a2ae208SSatish Balay #undef __FUNCT__ 12414a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 12428f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1243ff14e315SSatish Balay { 12443a40ed3dSBarry Smith PetscFunctionBegin; 124509b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12463a40ed3dSBarry Smith PetscFunctionReturn(0); 1247ff14e315SSatish Balay } 12480754003eSLois Curfman McInnes 12494a2ae208SSatish Balay #undef __FUNCT__ 12504a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 12517b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12520754003eSLois Curfman McInnes { 1253c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1254273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 1255182d2002SSatish Balay Scalar *av,*bv,*v = mat->v; 12560754003eSLois Curfman McInnes Mat newmat; 12570754003eSLois Curfman McInnes 12583a40ed3dSBarry Smith PetscFunctionBegin; 125978b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 126078b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1261e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1262e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 12630754003eSLois Curfman McInnes 1264182d2002SSatish Balay /* Check submatrixcall */ 1265182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1266182d2002SSatish Balay int n_cols,n_rows; 1267182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 126829bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1269182d2002SSatish Balay newmat = *B; 1270182d2002SSatish Balay } else { 12710754003eSLois Curfman McInnes /* Create and fill new matrix */ 1272b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1273182d2002SSatish Balay } 1274182d2002SSatish Balay 1275182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1276182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1277182d2002SSatish Balay 1278182d2002SSatish Balay for (i=0; i<ncols; i++) { 1279182d2002SSatish Balay av = v + m*icol[i]; 1280182d2002SSatish Balay for (j=0; j<nrows; j++) { 1281182d2002SSatish Balay *bv++ = av[irow[j]]; 12820754003eSLois Curfman McInnes } 12830754003eSLois Curfman McInnes } 1284182d2002SSatish Balay 1285182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 12866d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12876d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12880754003eSLois Curfman McInnes 12890754003eSLois Curfman McInnes /* Free work space */ 129078b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 129178b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1292182d2002SSatish Balay *B = newmat; 12933a40ed3dSBarry Smith PetscFunctionReturn(0); 12940754003eSLois Curfman McInnes } 12950754003eSLois Curfman McInnes 12964a2ae208SSatish Balay #undef __FUNCT__ 12974a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 12987b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1299905e6a2fSBarry Smith { 1300905e6a2fSBarry Smith int ierr,i; 1301905e6a2fSBarry Smith 13023a40ed3dSBarry Smith PetscFunctionBegin; 1303905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1304b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1305905e6a2fSBarry Smith } 1306905e6a2fSBarry Smith 1307905e6a2fSBarry Smith for (i=0; i<n; i++) { 13086a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1309905e6a2fSBarry Smith } 13103a40ed3dSBarry Smith PetscFunctionReturn(0); 1311905e6a2fSBarry Smith } 1312905e6a2fSBarry Smith 13134a2ae208SSatish Balay #undef __FUNCT__ 13144a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1315cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 13164b0e389bSBarry Smith { 13174b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 13183a40ed3dSBarry Smith int ierr; 1319273d9f13SBarry Smith PetscTruth flag; 13203a40ed3dSBarry Smith 13213a40ed3dSBarry Smith PetscFunctionBegin; 1322273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr); 1323273d9f13SBarry Smith if (!flag) { 1324cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 13253a40ed3dSBarry Smith PetscFunctionReturn(0); 13263a40ed3dSBarry Smith } 1327273d9f13SBarry Smith if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1328273d9f13SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(Scalar));CHKERRQ(ierr); 1329273d9f13SBarry Smith PetscFunctionReturn(0); 1330273d9f13SBarry Smith } 1331273d9f13SBarry Smith 13324a2ae208SSatish Balay #undef __FUNCT__ 13334a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1334273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1335273d9f13SBarry Smith { 1336273d9f13SBarry Smith int ierr; 1337273d9f13SBarry Smith 1338273d9f13SBarry Smith PetscFunctionBegin; 1339273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 13403a40ed3dSBarry Smith PetscFunctionReturn(0); 13414b0e389bSBarry Smith } 13424b0e389bSBarry Smith 1343289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1344a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1345905e6a2fSBarry Smith MatGetRow_SeqDense, 1346905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1347905e6a2fSBarry Smith MatMult_SeqDense, 1348905e6a2fSBarry Smith MatMultAdd_SeqDense, 13497c922b88SBarry Smith MatMultTranspose_SeqDense, 13507c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1351905e6a2fSBarry Smith MatSolve_SeqDense, 1352905e6a2fSBarry Smith MatSolveAdd_SeqDense, 13537c922b88SBarry Smith MatSolveTranspose_SeqDense, 13547c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1355905e6a2fSBarry Smith MatLUFactor_SeqDense, 1356905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1357ec8511deSBarry Smith MatRelax_SeqDense, 1358ec8511deSBarry Smith MatTranspose_SeqDense, 1359905e6a2fSBarry Smith MatGetInfo_SeqDense, 1360905e6a2fSBarry Smith MatEqual_SeqDense, 1361905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1362905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1363905e6a2fSBarry Smith MatNorm_SeqDense, 1364905e6a2fSBarry Smith 0, 1365905e6a2fSBarry Smith 0, 1366905e6a2fSBarry Smith 0, 1367905e6a2fSBarry Smith MatSetOption_SeqDense, 1368905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1369905e6a2fSBarry Smith MatZeroRows_SeqDense, 1370905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1371905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1372905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1373905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1374273d9f13SBarry Smith MatSetUpPreallocation_SeqDense, 1375273d9f13SBarry Smith 0, 1376905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1377905e6a2fSBarry Smith 0, 1378905e6a2fSBarry Smith 0, 1379905e6a2fSBarry Smith MatGetArray_SeqDense, 1380905e6a2fSBarry Smith MatRestoreArray_SeqDense, 13815609ef8eSBarry Smith MatDuplicate_SeqDense, 1382a5ae1ecdSBarry Smith 0, 1383a5ae1ecdSBarry Smith 0, 1384a5ae1ecdSBarry Smith 0, 1385a5ae1ecdSBarry Smith 0, 1386a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1387a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1388a5ae1ecdSBarry Smith 0, 13894b0e389bSBarry Smith MatGetValues_SeqDense, 1390a5ae1ecdSBarry Smith MatCopy_SeqDense, 1391a5ae1ecdSBarry Smith 0, 1392a5ae1ecdSBarry Smith MatScale_SeqDense, 1393a5ae1ecdSBarry Smith 0, 1394a5ae1ecdSBarry Smith 0, 1395a5ae1ecdSBarry Smith 0, 1396a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1397a5ae1ecdSBarry Smith 0, 1398a5ae1ecdSBarry Smith 0, 1399a5ae1ecdSBarry Smith 0, 1400a5ae1ecdSBarry Smith 0, 1401a5ae1ecdSBarry Smith 0, 1402a5ae1ecdSBarry Smith 0, 1403a5ae1ecdSBarry Smith 0, 1404a5ae1ecdSBarry Smith 0, 1405a5ae1ecdSBarry Smith 0, 1406a5ae1ecdSBarry Smith 0, 1407e03a110bSBarry Smith MatDestroy_SeqDense, 1408e03a110bSBarry Smith MatView_SeqDense, 1409a5ae1ecdSBarry Smith MatGetMaps_Petsc}; 141090ace30eSBarry Smith 14114a2ae208SSatish Balay #undef __FUNCT__ 14124a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 14134b828684SBarry Smith /*@C 1414fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1415d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1416d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1417289bc588SBarry Smith 1418db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1419db81eaa0SLois Curfman McInnes 142020563c6bSBarry Smith Input Parameters: 1421db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 14220c775827SLois Curfman McInnes . m - number of rows 142318f449edSLois Curfman McInnes . n - number of columns 1424db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1425dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 142620563c6bSBarry Smith 142720563c6bSBarry Smith Output Parameter: 142844cd7ae7SLois Curfman McInnes . A - the matrix 142920563c6bSBarry Smith 1430b259b22eSLois Curfman McInnes Notes: 143118f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 143218f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1433b4fd4287SBarry Smith set data=PETSC_NULL. 143418f449edSLois Curfman McInnes 1435027ccd11SLois Curfman McInnes Level: intermediate 1436027ccd11SLois Curfman McInnes 1437dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1438d65003e9SLois Curfman McInnes 1439db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 144020563c6bSBarry Smith @*/ 144144cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1442289bc588SBarry Smith { 1443273d9f13SBarry Smith int ierr; 14443b2fbd54SBarry Smith 14453a40ed3dSBarry Smith PetscFunctionBegin; 1446273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1447273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1448273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1449273d9f13SBarry Smith PetscFunctionReturn(0); 1450273d9f13SBarry Smith } 1451273d9f13SBarry Smith 14524a2ae208SSatish Balay #undef __FUNCT__ 14534a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1454273d9f13SBarry Smith /*@C 1455273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1456273d9f13SBarry Smith 1457273d9f13SBarry Smith Collective on MPI_Comm 1458273d9f13SBarry Smith 1459273d9f13SBarry Smith Input Parameters: 1460273d9f13SBarry Smith + A - the matrix 1461273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1462273d9f13SBarry Smith 1463273d9f13SBarry Smith Notes: 1464273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1465273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1466273d9f13SBarry Smith set data=PETSC_NULL. 1467273d9f13SBarry Smith 1468273d9f13SBarry Smith Level: intermediate 1469273d9f13SBarry Smith 1470273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1471273d9f13SBarry Smith 1472273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1473273d9f13SBarry Smith @*/ 1474273d9f13SBarry Smith int MatSeqDenseSetPreallocation(Mat B,Scalar *data) 1475273d9f13SBarry Smith { 1476273d9f13SBarry Smith Mat_SeqDense *b; 1477273d9f13SBarry Smith int ierr; 1478273d9f13SBarry Smith PetscTruth flg2; 1479273d9f13SBarry Smith 1480273d9f13SBarry Smith PetscFunctionBegin; 1481273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr); 1482273d9f13SBarry Smith if (!flg2) PetscFunctionReturn(0); 1483273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1484273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1485273d9f13SBarry Smith if (!data) { 1486b0a32e0cSBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(Scalar),&b->v);CHKERRQ(ierr); 1487273d9f13SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(Scalar));CHKERRQ(ierr); 1488273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 1489b0a32e0cSBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(Scalar)); 1490273d9f13SBarry Smith } else { /* user-allocated storage */ 1491273d9f13SBarry Smith b->v = data; 1492273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1493273d9f13SBarry Smith } 1494273d9f13SBarry Smith PetscFunctionReturn(0); 1495273d9f13SBarry Smith } 1496273d9f13SBarry Smith 1497273d9f13SBarry Smith EXTERN_C_BEGIN 14984a2ae208SSatish Balay #undef __FUNCT__ 14994a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1500273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1501273d9f13SBarry Smith { 1502273d9f13SBarry Smith Mat_SeqDense *b; 1503273d9f13SBarry Smith int ierr,size; 1504273d9f13SBarry Smith 1505273d9f13SBarry Smith PetscFunctionBegin; 1506273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 150729bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 150855659b69SBarry Smith 1509273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1510273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1511273d9f13SBarry Smith 1512b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1513549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1514549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 151544cd7ae7SLois Curfman McInnes B->factor = 0; 151690f02eecSBarry Smith B->mapping = 0; 1517b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 151844cd7ae7SLois Curfman McInnes B->data = (void*)b; 151918f449edSLois Curfman McInnes 1520273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 1521273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1522a5ae1ecdSBarry Smith 152344cd7ae7SLois Curfman McInnes b->pivots = 0; 1524273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1525273d9f13SBarry Smith b->v = 0; 15264e220ebcSLois Curfman McInnes 15273a40ed3dSBarry Smith PetscFunctionReturn(0); 1528289bc588SBarry Smith } 1529273d9f13SBarry Smith EXTERN_C_END 15303b2fbd54SBarry Smith 15313b2fbd54SBarry Smith 15323b2fbd54SBarry Smith 15333b2fbd54SBarry Smith 1534