173f4d377SMatthew Knepley /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/ 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" 11607cd303SBarry Smith int MatAXPY_SeqDense(PetscScalar *alpha,Mat X,Mat Y,MatStructure str) 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; 2887828ca2SBarry Smith PetscScalar *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" 5387828ca2SBarry Smith int MatScale_SeqDense(PetscScalar *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); 82ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRF) 83ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable."); 84ae7cfcebSSatish Balay #else 85273d9f13SBarry Smith LAgetrf_(&A->m,&A->n,mat->v,&A->m,mat->pivots,&info); 8629bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 8729bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 882ffef20aSBarry Smith #endif 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) { 10587828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));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 */ 13287828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));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" 140*15e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,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__ "MatCholeskyFactor_SeqDense" 151*15e8a5b3SHong Zhang int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 152289bc588SBarry Smith { 153c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 154606d414cSSatish Balay int info,ierr; 15555659b69SBarry Smith 1563a40ed3dSBarry Smith PetscFunctionBegin; 157752f5567SLois Curfman McInnes if (mat->pivots) { 158606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 159b0a32e0cSBarry Smith PetscLogObjectMemory(A,-A->m*sizeof(int)); 160752f5567SLois Curfman McInnes mat->pivots = 0; 161752f5567SLois Curfman McInnes } 162273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 163ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRF) 164ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable."); 165ae7cfcebSSatish Balay #else 166273d9f13SBarry Smith LApotrf_("L",&A->n,mat->v,&A->m,&info); 16729bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 1682ffef20aSBarry Smith #endif 169c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 170b0a32e0cSBarry Smith PetscLogFlops((A->n*A->n*A->n)/3); 1713a40ed3dSBarry Smith PetscFunctionReturn(0); 172289bc588SBarry Smith } 173289bc588SBarry Smith 1744a2ae208SSatish Balay #undef __FUNCT__ 1750b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 1760b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 1770b4b3355SBarry Smith { 1780b4b3355SBarry Smith int ierr; 179*15e8a5b3SHong Zhang MatFactorInfo info; 1800b4b3355SBarry Smith 1810b4b3355SBarry Smith PetscFunctionBegin; 182*15e8a5b3SHong Zhang info.fill = 1.0; 183*15e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 1840b4b3355SBarry Smith PetscFunctionReturn(0); 1850b4b3355SBarry Smith } 1860b4b3355SBarry Smith 1870b4b3355SBarry Smith #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 1898f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 190289bc588SBarry Smith { 191c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 192a57ad546SLois Curfman McInnes int one = 1,info,ierr; 19387828ca2SBarry Smith PetscScalar *x,*y; 19467e560aaSBarry Smith 1953a40ed3dSBarry Smith PetscFunctionBegin; 196273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 1977a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1987a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 19987828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 200c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 201ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 202ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 203ae7cfcebSSatish Balay #else 204273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 2052ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 206ae7cfcebSSatish Balay #endif 2077a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 208ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 209ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 210ae7cfcebSSatish Balay #else 211273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 2122ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 213ae7cfcebSSatish Balay #endif 214289bc588SBarry Smith } 21529bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2167a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2177a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 218b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2193a40ed3dSBarry Smith PetscFunctionReturn(0); 220289bc588SBarry Smith } 2216ee01492SSatish Balay 2224a2ae208SSatish Balay #undef __FUNCT__ 2234a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 2247c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 225da3a660dSBarry Smith { 226c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2277a97a34bSBarry Smith int ierr,one = 1,info; 22887828ca2SBarry Smith PetscScalar *x,*y; 22967e560aaSBarry Smith 2303a40ed3dSBarry Smith PetscFunctionBegin; 231273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2327a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2337a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 23487828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 235752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 236da3a660dSBarry Smith if (mat->pivots) { 237ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 238ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 239ae7cfcebSSatish Balay #else 240273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 2412ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 242ae7cfcebSSatish Balay #endif 2437a97a34bSBarry Smith } else { 244ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 245ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 246ae7cfcebSSatish Balay #else 247273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 2482ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 249ae7cfcebSSatish Balay #endif 250da3a660dSBarry Smith } 2517a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2527a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 253b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2543a40ed3dSBarry Smith PetscFunctionReturn(0); 255da3a660dSBarry Smith } 2566ee01492SSatish Balay 2574a2ae208SSatish Balay #undef __FUNCT__ 2584a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 2598f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 260da3a660dSBarry Smith { 261c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2627a97a34bSBarry Smith int ierr,one = 1,info; 26387828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 264da3a660dSBarry Smith Vec tmp = 0; 26567e560aaSBarry Smith 2663a40ed3dSBarry Smith PetscFunctionBegin; 2677a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2687a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 269273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 270da3a660dSBarry Smith if (yy == zz) { 27178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 272b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 27378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 274da3a660dSBarry Smith } 27587828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 276752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 277da3a660dSBarry Smith if (mat->pivots) { 278ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 279ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 280ae7cfcebSSatish Balay #else 281273d9f13SBarry Smith LAgetrs_("N",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 2822ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 283ae7cfcebSSatish Balay #endif 284a8c6a408SBarry Smith } else { 285ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 286ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 287ae7cfcebSSatish Balay #else 288273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 2892ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 290ae7cfcebSSatish Balay #endif 291da3a660dSBarry Smith } 292a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 293a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 2947a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2957a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 296b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 298da3a660dSBarry Smith } 29967e560aaSBarry Smith 3004a2ae208SSatish Balay #undef __FUNCT__ 3014a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 3027c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 303da3a660dSBarry Smith { 304c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3056abc6512SBarry Smith int one = 1,info,ierr; 30687828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 307da3a660dSBarry Smith Vec tmp; 30867e560aaSBarry Smith 3093a40ed3dSBarry Smith PetscFunctionBegin; 310273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3117a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3127a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 313da3a660dSBarry Smith if (yy == zz) { 31478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 315b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 31678b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 317da3a660dSBarry Smith } 31887828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 319752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 320da3a660dSBarry Smith if (mat->pivots) { 321ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 322ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 323ae7cfcebSSatish Balay #else 324273d9f13SBarry Smith LAgetrs_("T",&A->m,&one,mat->v,&A->m,mat->pivots,y,&A->m,&info); 3252ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 326ae7cfcebSSatish Balay #endif 3273a40ed3dSBarry Smith } else { 328ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 329ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 330ae7cfcebSSatish Balay #else 331273d9f13SBarry Smith LApotrs_("L",&A->m,&one,mat->v,&A->m,y,&A->m,&info); 3322ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 333ae7cfcebSSatish Balay #endif 334da3a660dSBarry Smith } 33590f02eecSBarry Smith if (tmp) { 33690f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 33790f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3383a40ed3dSBarry Smith } else { 33990f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 34090f02eecSBarry Smith } 3417a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3427a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 343b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3443a40ed3dSBarry Smith PetscFunctionReturn(0); 345da3a660dSBarry Smith } 346289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3474a2ae208SSatish Balay #undef __FUNCT__ 3484a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 349329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 350c14dc6b6SHong Zhang PetscReal shift,int its,int lits,Vec xx) 351289bc588SBarry Smith { 352c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 35387828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 354273d9f13SBarry Smith int ierr,m = A->m,i; 355aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 356bc1b551cSSatish Balay int o = 1; 357bc1b551cSSatish Balay #endif 358289bc588SBarry Smith 3593a40ed3dSBarry Smith PetscFunctionBegin; 360289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3613a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 362bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 363289bc588SBarry Smith } 3647a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3657a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 366b965ef7fSBarry Smith its = its*lits; 36791723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 368289bc588SBarry Smith while (its--) { 369289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 370289bc588SBarry Smith for (i=0; i<m; i++) { 371aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 372f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 373f1747703SBarry Smith not happy about returning a double complex */ 374f1747703SBarry Smith int _i; 37587828ca2SBarry Smith PetscScalar sum = b[i]; 376f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3773f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 378f1747703SBarry Smith } 379f1747703SBarry Smith xt = sum; 380f1747703SBarry Smith #else 381289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 382f1747703SBarry Smith #endif 38355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 384289bc588SBarry Smith } 385289bc588SBarry Smith } 386289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 387289bc588SBarry Smith for (i=m-1; i>=0; i--) { 388aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 389f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 390f1747703SBarry Smith not happy about returning a double complex */ 391f1747703SBarry Smith int _i; 39287828ca2SBarry Smith PetscScalar sum = b[i]; 393f1747703SBarry Smith for (_i=0; _i<m; _i++) { 3943f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 395f1747703SBarry Smith } 396f1747703SBarry Smith xt = sum; 397f1747703SBarry Smith #else 398289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 399f1747703SBarry Smith #endif 40055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 401289bc588SBarry Smith } 402289bc588SBarry Smith } 403289bc588SBarry Smith } 404600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4057a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4063a40ed3dSBarry Smith PetscFunctionReturn(0); 407289bc588SBarry Smith } 408289bc588SBarry Smith 409289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4104a2ae208SSatish Balay #undef __FUNCT__ 4114a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 4127c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 413289bc588SBarry Smith { 414c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 41587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 416ea709b57SSatish Balay int ierr,_One=1; 417ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4183a40ed3dSBarry Smith 4193a40ed3dSBarry Smith PetscFunctionBegin; 420273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4217a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4227a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 423273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 4247a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4257a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 426b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->n); 4273a40ed3dSBarry Smith PetscFunctionReturn(0); 428289bc588SBarry Smith } 4296ee01492SSatish Balay 4304a2ae208SSatish Balay #undef __FUNCT__ 4314a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 43244cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 433289bc588SBarry Smith { 434c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 43587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 436329f5518SBarry Smith int ierr,_One=1; 4373a40ed3dSBarry Smith 4383a40ed3dSBarry Smith PetscFunctionBegin; 439273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4407a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4417a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 442273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DZero,y,&_One); 4437a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4447a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 445b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->m); 4463a40ed3dSBarry Smith PetscFunctionReturn(0); 447289bc588SBarry Smith } 4486ee01492SSatish Balay 4494a2ae208SSatish Balay #undef __FUNCT__ 4504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 45144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 452289bc588SBarry Smith { 453c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 45487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 455329f5518SBarry Smith int ierr,_One=1; 4563a40ed3dSBarry Smith 4573a40ed3dSBarry Smith PetscFunctionBegin; 458273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 459600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4607a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4617a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 462273d9f13SBarry Smith LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4637a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4647a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 465b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4663a40ed3dSBarry Smith PetscFunctionReturn(0); 467289bc588SBarry Smith } 4686ee01492SSatish Balay 4694a2ae208SSatish Balay #undef __FUNCT__ 4704a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 4717c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 472289bc588SBarry Smith { 473c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 47487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 4757a97a34bSBarry Smith int ierr,_One=1; 47687828ca2SBarry Smith PetscScalar _DOne=1.0; 4773a40ed3dSBarry Smith 4783a40ed3dSBarry Smith PetscFunctionBegin; 479273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 480600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4817a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4827a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 483273d9f13SBarry Smith LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(A->m),x,&_One,&_DOne,y,&_One); 4847a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4857a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 486b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4873a40ed3dSBarry Smith PetscFunctionReturn(0); 488289bc588SBarry Smith } 489289bc588SBarry Smith 490289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4914a2ae208SSatish Balay #undef __FUNCT__ 4924a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 49387828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 494289bc588SBarry Smith { 495c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 49687828ca2SBarry Smith PetscScalar *v; 497b0a32e0cSBarry Smith int i,ierr; 49867e560aaSBarry Smith 4993a40ed3dSBarry Smith PetscFunctionBegin; 500273d9f13SBarry Smith *ncols = A->n; 501289bc588SBarry Smith if (cols) { 502b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 503273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 504289bc588SBarry Smith } 505289bc588SBarry Smith if (vals) { 50687828ca2SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 507289bc588SBarry Smith v = mat->v + row; 508273d9f13SBarry Smith for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += A->m;} 509289bc588SBarry Smith } 5103a40ed3dSBarry Smith PetscFunctionReturn(0); 511289bc588SBarry Smith } 5126ee01492SSatish Balay 5134a2ae208SSatish Balay #undef __FUNCT__ 5144a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 51587828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 516289bc588SBarry Smith { 517606d414cSSatish Balay int ierr; 518606d414cSSatish Balay PetscFunctionBegin; 519606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 520606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5213a40ed3dSBarry Smith PetscFunctionReturn(0); 522289bc588SBarry Smith } 523289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5244a2ae208SSatish Balay #undef __FUNCT__ 5254a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 5268f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 52787828ca2SBarry Smith int *indexn,PetscScalar *v,InsertMode addv) 528289bc588SBarry Smith { 529c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 530289bc588SBarry Smith int i,j; 531d6dfbf8fSBarry Smith 5323a40ed3dSBarry Smith PetscFunctionBegin; 533289bc588SBarry Smith if (!mat->roworiented) { 534dbb450caSBarry Smith if (addv == INSERT_VALUES) { 535289bc588SBarry Smith for (j=0; j<n; j++) { 5365ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; 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 540289bc588SBarry Smith for (i=0; i<m; i++) { 5415ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 542aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 54329bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 54458804f6dSBarry Smith #endif 545273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] = *v++; 546289bc588SBarry Smith } 547289bc588SBarry Smith } 5483a40ed3dSBarry Smith } else { 549289bc588SBarry Smith for (j=0; j<n; j++) { 5505ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 551aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 55229bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 55358804f6dSBarry Smith #endif 554289bc588SBarry Smith for (i=0; i<m; i++) { 5555ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 556aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 55729bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 55858804f6dSBarry Smith #endif 559273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] += *v++; 560289bc588SBarry Smith } 561289bc588SBarry Smith } 562289bc588SBarry Smith } 5633a40ed3dSBarry Smith } else { 564dbb450caSBarry Smith if (addv == INSERT_VALUES) { 565e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 5665ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 567aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 56829bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56958804f6dSBarry Smith #endif 570e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 5715ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 572aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57329bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 57458804f6dSBarry Smith #endif 575273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] = *v++; 576e8d4e0b9SBarry Smith } 577e8d4e0b9SBarry Smith } 5783a40ed3dSBarry Smith } else { 579289bc588SBarry Smith for (i=0; i<m; i++) { 5805ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 581aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 58229bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 58358804f6dSBarry Smith #endif 584289bc588SBarry Smith for (j=0; j<n; j++) { 5855ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 586aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 58729bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 58858804f6dSBarry Smith #endif 589273d9f13SBarry Smith mat->v[indexn[j]*A->m + indexm[i]] += *v++; 590289bc588SBarry Smith } 591289bc588SBarry Smith } 592289bc588SBarry Smith } 593e8d4e0b9SBarry Smith } 5943a40ed3dSBarry Smith PetscFunctionReturn(0); 595289bc588SBarry Smith } 596e8d4e0b9SBarry Smith 5974a2ae208SSatish Balay #undef __FUNCT__ 5984a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 59987828ca2SBarry Smith int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v) 600ae80bb75SLois Curfman McInnes { 601ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 602ae80bb75SLois Curfman McInnes int i,j; 60387828ca2SBarry Smith PetscScalar *vpt = v; 604ae80bb75SLois Curfman McInnes 6053a40ed3dSBarry Smith PetscFunctionBegin; 606ae80bb75SLois Curfman McInnes /* row-oriented output */ 607ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 608ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 609273d9f13SBarry Smith *vpt++ = mat->v[indexn[j]*A->m + indexm[i]]; 610ae80bb75SLois Curfman McInnes } 611ae80bb75SLois Curfman McInnes } 6123a40ed3dSBarry Smith PetscFunctionReturn(0); 613ae80bb75SLois Curfman McInnes } 614ae80bb75SLois Curfman McInnes 615289bc588SBarry Smith /* -----------------------------------------------------------------*/ 616289bc588SBarry Smith 617e090d566SSatish Balay #include "petscsys.h" 61888e32edaSLois Curfman McInnes 619273d9f13SBarry Smith EXTERN_C_BEGIN 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 622b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A) 62388e32edaSLois Curfman McInnes { 62488e32edaSLois Curfman McInnes Mat_SeqDense *a; 62588e32edaSLois Curfman McInnes Mat B; 62688e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 62788e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 62887828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 62919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 63088e32edaSLois Curfman McInnes 6313a40ed3dSBarry Smith PetscFunctionBegin; 632d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 63329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 634b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6350752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 636552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 63788e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 63888e32edaSLois Curfman McInnes 63990ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 64090ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 64190ace30eSBarry Smith B = *A; 64290ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6434a41a4d3SSatish Balay v = a->v; 6444a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6454a41a4d3SSatish Balay from row major to column major */ 64687828ca2SBarry Smith ierr = PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 64790ace30eSBarry Smith /* read in nonzero values */ 6484a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6494a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6504a41a4d3SSatish Balay for (j=0; j<N; j++) { 6514a41a4d3SSatish Balay for (i=0; i<M; i++) { 6524a41a4d3SSatish Balay *v++ =w[i*N+j]; 6534a41a4d3SSatish Balay } 6544a41a4d3SSatish Balay } 655606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6566d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6576d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 65890ace30eSBarry Smith } else { 65988e32edaSLois Curfman McInnes /* read row lengths */ 660b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 6610752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 66288e32edaSLois Curfman McInnes 66388e32edaSLois Curfman McInnes /* create our matrix */ 664b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 66588e32edaSLois Curfman McInnes B = *A; 66688e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 66788e32edaSLois Curfman McInnes v = a->v; 66888e32edaSLois Curfman McInnes 66988e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 670b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 671b0a32e0cSBarry Smith cols = scols; 6720752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 67387828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 674b0a32e0cSBarry Smith vals = svals; 6750752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 67688e32edaSLois Curfman McInnes 67788e32edaSLois Curfman McInnes /* insert into matrix */ 67888e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 67988e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 68088e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 68188e32edaSLois Curfman McInnes } 682606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 683606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 684606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 68588e32edaSLois Curfman McInnes 6866d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6876d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68890ace30eSBarry Smith } 6893a40ed3dSBarry Smith PetscFunctionReturn(0); 69088e32edaSLois Curfman McInnes } 691273d9f13SBarry Smith EXTERN_C_END 69288e32edaSLois Curfman McInnes 693e090d566SSatish Balay #include "petscsys.h" 694932b0c3eSLois Curfman McInnes 6954a2ae208SSatish Balay #undef __FUNCT__ 6964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 697b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 698289bc588SBarry Smith { 699932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 700fb9695e5SSatish Balay int ierr,i,j; 701fb9695e5SSatish Balay char *name; 70287828ca2SBarry Smith PetscScalar *v; 703f3ef73ceSBarry Smith PetscViewerFormat format; 704932b0c3eSLois Curfman McInnes 7053a40ed3dSBarry Smith PetscFunctionBegin; 706435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 707b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 708456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7093a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 710fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 711b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 712273d9f13SBarry Smith for (i=0; i<A->m; i++) { 71344cd7ae7SLois Curfman McInnes v = a->v + i; 714b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 715273d9f13SBarry Smith for (j=0; j<A->n; j++) { 716aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 717329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 71862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 719329f5518SBarry Smith } else if (PetscRealPart(*v)) { 72062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7216831982aSBarry Smith } 72280cd9d93SLois Curfman McInnes #else 7236831982aSBarry Smith if (*v) { 72462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 7256831982aSBarry Smith } 72680cd9d93SLois Curfman McInnes #endif 727273d9f13SBarry Smith v += A->m; 72880cd9d93SLois Curfman McInnes } 729b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 73080cd9d93SLois Curfman McInnes } 731b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7323a40ed3dSBarry Smith } else { 733b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 734aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 735ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 73647989497SBarry Smith /* determine if matrix has all real values */ 73747989497SBarry Smith v = a->v; 738273d9f13SBarry Smith for (i=0; i<A->m*A->n; i++) { 739ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 74047989497SBarry Smith } 74147989497SBarry Smith #endif 742fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7433a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 744b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 745fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 746fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 747ffac6cdbSBarry Smith } 748ffac6cdbSBarry Smith 749273d9f13SBarry Smith for (i=0; i<A->m; i++) { 750932b0c3eSLois Curfman McInnes v = a->v + i; 751273d9f13SBarry Smith for (j=0; j<A->n; j++) { 752aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 75347989497SBarry Smith if (allreal) { 754b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); v += A->m; 75547989497SBarry Smith } else { 756b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); v += A->m; 75747989497SBarry Smith } 758289bc588SBarry Smith #else 759b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += A->m; 760289bc588SBarry Smith #endif 761289bc588SBarry Smith } 762b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 763289bc588SBarry Smith } 764fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 765b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 766ffac6cdbSBarry Smith } 767b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 768da3a660dSBarry Smith } 769b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7703a40ed3dSBarry Smith PetscFunctionReturn(0); 771289bc588SBarry Smith } 772289bc588SBarry Smith 7734a2ae208SSatish Balay #undef __FUNCT__ 7744a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 775b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 776932b0c3eSLois Curfman McInnes { 777932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 778273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 77987828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 780f3ef73ceSBarry Smith PetscViewerFormat format; 781932b0c3eSLois Curfman McInnes 7823a40ed3dSBarry Smith PetscFunctionBegin; 783b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 78490ace30eSBarry Smith 785b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 786fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 78790ace30eSBarry Smith /* store the matrix as a dense matrix */ 788b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 789552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 79090ace30eSBarry Smith col_lens[1] = m; 79190ace30eSBarry Smith col_lens[2] = n; 79290ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 7930752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 794606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 79590ace30eSBarry Smith 79690ace30eSBarry Smith /* write out matrix, by rows */ 79787828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 79890ace30eSBarry Smith v = a->v; 79990ace30eSBarry Smith for (i=0; i<m; i++) { 80090ace30eSBarry Smith for (j=0; j<n; j++) { 80190ace30eSBarry Smith vals[i + j*m] = *v++; 80290ace30eSBarry Smith } 80390ace30eSBarry Smith } 8040752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 805606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 80690ace30eSBarry Smith } else { 807b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 808552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 809932b0c3eSLois Curfman McInnes col_lens[1] = m; 810932b0c3eSLois Curfman McInnes col_lens[2] = n; 811932b0c3eSLois Curfman McInnes col_lens[3] = nz; 812932b0c3eSLois Curfman McInnes 813932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 814932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8150752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 816932b0c3eSLois Curfman McInnes 817932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 818932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 819932b0c3eSLois Curfman McInnes ict = 0; 820932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 821932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 822932b0c3eSLois Curfman McInnes } 8230752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 824606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 825932b0c3eSLois Curfman McInnes 826932b0c3eSLois Curfman McInnes /* store nonzero values */ 82787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 828932b0c3eSLois Curfman McInnes ict = 0; 829932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 830932b0c3eSLois Curfman McInnes v = a->v + i; 831932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 832273d9f13SBarry Smith anonz[ict++] = *v; v += A->m; 833932b0c3eSLois Curfman McInnes } 834932b0c3eSLois Curfman McInnes } 8350752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 836606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 83790ace30eSBarry Smith } 8383a40ed3dSBarry Smith PetscFunctionReturn(0); 839932b0c3eSLois Curfman McInnes } 840932b0c3eSLois Curfman McInnes 8414a2ae208SSatish Balay #undef __FUNCT__ 8424a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 843b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 844f1af5d2fSBarry Smith { 845f1af5d2fSBarry Smith Mat A = (Mat) Aa; 846f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 847fb9695e5SSatish Balay int m = A->m,n = A->n,color,i,j,ierr; 84887828ca2SBarry Smith PetscScalar *v = a->v; 849b0a32e0cSBarry Smith PetscViewer viewer; 850b0a32e0cSBarry Smith PetscDraw popup; 851329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 852f3ef73ceSBarry Smith PetscViewerFormat format; 853f1af5d2fSBarry Smith 854f1af5d2fSBarry Smith PetscFunctionBegin; 855f1af5d2fSBarry Smith 856f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 857b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 858b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 859f1af5d2fSBarry Smith 860f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 861fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 862f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 863b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 864f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 865f1af5d2fSBarry Smith x_l = j; 866f1af5d2fSBarry Smith x_r = x_l + 1.0; 867f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 868f1af5d2fSBarry Smith y_l = m - i - 1.0; 869f1af5d2fSBarry Smith y_r = y_l + 1.0; 870f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 871329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 872b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 873329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 874b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 875f1af5d2fSBarry Smith } else { 876f1af5d2fSBarry Smith continue; 877f1af5d2fSBarry Smith } 878f1af5d2fSBarry Smith #else 879f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 880b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 881f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 882b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 883f1af5d2fSBarry Smith } else { 884f1af5d2fSBarry Smith continue; 885f1af5d2fSBarry Smith } 886f1af5d2fSBarry Smith #endif 887b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 888f1af5d2fSBarry Smith } 889f1af5d2fSBarry Smith } 890f1af5d2fSBarry Smith } else { 891f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 892f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 893f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 894f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 895f1af5d2fSBarry Smith } 896b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 897b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 898b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 899f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 900f1af5d2fSBarry Smith x_l = j; 901f1af5d2fSBarry Smith x_r = x_l + 1.0; 902f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 903f1af5d2fSBarry Smith y_l = m - i - 1.0; 904f1af5d2fSBarry Smith y_r = y_l + 1.0; 905b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 906b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 907f1af5d2fSBarry Smith } 908f1af5d2fSBarry Smith } 909f1af5d2fSBarry Smith } 910f1af5d2fSBarry Smith PetscFunctionReturn(0); 911f1af5d2fSBarry Smith } 912f1af5d2fSBarry Smith 9134a2ae208SSatish Balay #undef __FUNCT__ 9144a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 915b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 916f1af5d2fSBarry Smith { 917b0a32e0cSBarry Smith PetscDraw draw; 918f1af5d2fSBarry Smith PetscTruth isnull; 919329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 920f1af5d2fSBarry Smith int ierr; 921f1af5d2fSBarry Smith 922f1af5d2fSBarry Smith PetscFunctionBegin; 923b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 924b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 925f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 926f1af5d2fSBarry Smith 927f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 928273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 929f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 930b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 931b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 932f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 933f1af5d2fSBarry Smith PetscFunctionReturn(0); 934f1af5d2fSBarry Smith } 935f1af5d2fSBarry Smith 9364a2ae208SSatish Balay #undef __FUNCT__ 9374a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 938b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer) 939932b0c3eSLois Curfman McInnes { 940932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 941bcd2baecSBarry Smith int ierr; 942f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 943932b0c3eSLois Curfman McInnes 9443a40ed3dSBarry Smith PetscFunctionBegin; 945b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 946b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 947fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 948fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 9490f5bd95cSBarry Smith 9500f5bd95cSBarry Smith if (issocket) { 951b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9520f5bd95cSBarry Smith } else if (isascii) { 9533a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9540f5bd95cSBarry Smith } else if (isbinary) { 9553a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 956f1af5d2fSBarry Smith } else if (isdraw) { 957f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9585cd90555SBarry Smith } else { 95929bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 960932b0c3eSLois Curfman McInnes } 9613a40ed3dSBarry Smith PetscFunctionReturn(0); 962932b0c3eSLois Curfman McInnes } 963289bc588SBarry Smith 9644a2ae208SSatish Balay #undef __FUNCT__ 9654a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 966c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 967289bc588SBarry Smith { 968ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 96990f02eecSBarry Smith int ierr; 97090f02eecSBarry Smith 9713a40ed3dSBarry Smith PetscFunctionBegin; 972aa482453SBarry Smith #if defined(PETSC_USE_LOG) 973b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 974a5a9c739SBarry Smith #endif 975606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 976606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 977606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 9783a40ed3dSBarry Smith PetscFunctionReturn(0); 979289bc588SBarry Smith } 980289bc588SBarry Smith 9814a2ae208SSatish Balay #undef __FUNCT__ 9824a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 9838f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 984289bc588SBarry Smith { 985c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 986549d3d68SSatish Balay int k,j,m,n,ierr; 98787828ca2SBarry Smith PetscScalar *v,tmp; 98848b35521SBarry Smith 9893a40ed3dSBarry Smith PetscFunctionBegin; 990273d9f13SBarry Smith v = mat->v; m = A->m; n = A->n; 9917c922b88SBarry Smith if (!matout) { /* in place transpose */ 992d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 99387828ca2SBarry Smith PetscScalar *w; 99487828ca2SBarry Smith ierr = PetscMalloc((m+1)*(n+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 995d64ed03dSBarry Smith for (j=0; j<m; j++) { 996d64ed03dSBarry Smith for (k=0; k<n; k++) { 997d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 998d64ed03dSBarry Smith } 999d64ed03dSBarry Smith } 100087828ca2SBarry Smith ierr = PetscMemcpy(v,w,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1001606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 1002d64ed03dSBarry Smith } else { 1003d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1004289bc588SBarry Smith for (k=0; k<j; k++) { 1005d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 1006d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 1007d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 1008289bc588SBarry Smith } 1009289bc588SBarry Smith } 1010d64ed03dSBarry Smith } 10113a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1012d3e5ee88SLois Curfman McInnes Mat tmat; 1013ec8511deSBarry Smith Mat_SeqDense *tmatd; 101487828ca2SBarry Smith PetscScalar *v2; 1015ea709b57SSatish Balay 1016273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1017ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10180de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1019d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10200de55854SLois Curfman McInnes for (k=0; k<m; k++) v2[j + k*n] = v[k + j*m]; 1021d3e5ee88SLois Curfman McInnes } 10226d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10236d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1024d3e5ee88SLois Curfman McInnes *matout = tmat; 102548b35521SBarry Smith } 10263a40ed3dSBarry Smith PetscFunctionReturn(0); 1027289bc588SBarry Smith } 1028289bc588SBarry Smith 10294a2ae208SSatish Balay #undef __FUNCT__ 10304a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 10318f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1032289bc588SBarry Smith { 1033c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1034c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1035273d9f13SBarry Smith int ierr,i; 103687828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1037273d9f13SBarry Smith PetscTruth flag; 10389ea5d5aeSSatish Balay 10393a40ed3dSBarry Smith PetscFunctionBegin; 1040273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr); 1041273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 1042273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1043273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1044273d9f13SBarry Smith for (i=0; i<A1->m*A1->n; i++) { 10453a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1046289bc588SBarry Smith v1++; v2++; 1047289bc588SBarry Smith } 104877c4ece6SBarry Smith *flg = PETSC_TRUE; 10493a40ed3dSBarry Smith PetscFunctionReturn(0); 1050289bc588SBarry Smith } 1051289bc588SBarry Smith 10524a2ae208SSatish Balay #undef __FUNCT__ 10534a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 10548f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1055289bc588SBarry Smith { 1056c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10577a97a34bSBarry Smith int ierr,i,n,len; 105887828ca2SBarry Smith PetscScalar *x,zero = 0.0; 105944cd7ae7SLois Curfman McInnes 10603a40ed3dSBarry Smith PetscFunctionBegin; 10617a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10627a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1063600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1064273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1065273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 106644cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 1067273d9f13SBarry Smith x[i] = mat->v[i*A->m + i]; 1068289bc588SBarry Smith } 10697a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10703a40ed3dSBarry Smith PetscFunctionReturn(0); 1071289bc588SBarry Smith } 1072289bc588SBarry Smith 10734a2ae208SSatish Balay #undef __FUNCT__ 10744a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 10758f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1076289bc588SBarry Smith { 1077c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 107887828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1079273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 108055659b69SBarry Smith 10813a40ed3dSBarry Smith PetscFunctionBegin; 108228988994SBarry Smith if (ll) { 10837a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1084600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1085273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1086da3a660dSBarry Smith for (i=0; i<m; i++) { 1087da3a660dSBarry Smith x = l[i]; 1088da3a660dSBarry Smith v = mat->v + i; 1089da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1090da3a660dSBarry Smith } 10917a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1092b0a32e0cSBarry Smith PetscLogFlops(n*m); 1093da3a660dSBarry Smith } 109428988994SBarry Smith if (rr) { 10957a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1096600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1097273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1098da3a660dSBarry Smith for (i=0; i<n; i++) { 1099da3a660dSBarry Smith x = r[i]; 1100da3a660dSBarry Smith v = mat->v + i*m; 1101da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1102da3a660dSBarry Smith } 11037a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1104b0a32e0cSBarry Smith PetscLogFlops(n*m); 1105da3a660dSBarry Smith } 11063a40ed3dSBarry Smith PetscFunctionReturn(0); 1107289bc588SBarry Smith } 1108289bc588SBarry Smith 11094a2ae208SSatish Balay #undef __FUNCT__ 11104a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1111064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1112289bc588SBarry Smith { 1113c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 111487828ca2SBarry Smith PetscScalar *v = mat->v; 1115329f5518SBarry Smith PetscReal sum = 0.0; 1116289bc588SBarry Smith int i,j; 111755659b69SBarry Smith 11183a40ed3dSBarry Smith PetscFunctionBegin; 1119289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1120273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1121aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1122329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1123289bc588SBarry Smith #else 1124289bc588SBarry Smith sum += (*v)*(*v); v++; 1125289bc588SBarry Smith #endif 1126289bc588SBarry Smith } 1127064f8208SBarry Smith *nrm = sqrt(sum); 1128b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11293a40ed3dSBarry Smith } else if (type == NORM_1) { 1130064f8208SBarry Smith *nrm = 0.0; 1131273d9f13SBarry Smith for (j=0; j<A->n; j++) { 1132289bc588SBarry Smith sum = 0.0; 1133273d9f13SBarry Smith for (i=0; i<A->m; i++) { 113433a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1135289bc588SBarry Smith } 1136064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1137289bc588SBarry Smith } 1138b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11393a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1140064f8208SBarry Smith *nrm = 0.0; 1141273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1142289bc588SBarry Smith v = mat->v + j; 1143289bc588SBarry Smith sum = 0.0; 1144273d9f13SBarry Smith for (i=0; i<A->n; i++) { 1145273d9f13SBarry Smith sum += PetscAbsScalar(*v); v += A->m; 1146289bc588SBarry Smith } 1147064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1148289bc588SBarry Smith } 1149b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11503a40ed3dSBarry Smith } else { 115129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1152289bc588SBarry Smith } 11533a40ed3dSBarry Smith PetscFunctionReturn(0); 1154289bc588SBarry Smith } 1155289bc588SBarry Smith 11564a2ae208SSatish Balay #undef __FUNCT__ 11574a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 11588f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1159289bc588SBarry Smith { 1160c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 116167e560aaSBarry Smith 11623a40ed3dSBarry Smith PetscFunctionBegin; 1163b5a2b587SKris Buschelman switch (op) { 1164b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1165b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1166b5a2b587SKris Buschelman break; 1167b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1168b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1169b5a2b587SKris Buschelman break; 1170b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1171b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1172b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1173b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1174b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1175b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1176b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1177b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1178b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1179b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1180b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1181b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1182b5a2b587SKris Buschelman break; 1183b5a2b587SKris Buschelman default: 118429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 11853a40ed3dSBarry Smith } 11863a40ed3dSBarry Smith PetscFunctionReturn(0); 1187289bc588SBarry Smith } 1188289bc588SBarry Smith 11894a2ae208SSatish Balay #undef __FUNCT__ 11904a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 11918f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 11926f0a148fSBarry Smith { 1193ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1194549d3d68SSatish Balay int ierr; 11953a40ed3dSBarry Smith 11963a40ed3dSBarry Smith PetscFunctionBegin; 119787828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 11983a40ed3dSBarry Smith PetscFunctionReturn(0); 11996f0a148fSBarry Smith } 12006f0a148fSBarry Smith 12014a2ae208SSatish Balay #undef __FUNCT__ 12024a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 12038f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 12044e220ebcSLois Curfman McInnes { 12053a40ed3dSBarry Smith PetscFunctionBegin; 12064e220ebcSLois Curfman McInnes *bs = 1; 12073a40ed3dSBarry Smith PetscFunctionReturn(0); 12084e220ebcSLois Curfman McInnes } 12094e220ebcSLois Curfman McInnes 12104a2ae208SSatish Balay #undef __FUNCT__ 12114a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 121287828ca2SBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag) 12136f0a148fSBarry Smith { 1214ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1215273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 121687828ca2SBarry Smith PetscScalar *slot; 121755659b69SBarry Smith 12183a40ed3dSBarry Smith PetscFunctionBegin; 1219e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 122078b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12216f0a148fSBarry Smith for (i=0; i<N; i++) { 12226f0a148fSBarry Smith slot = l->v + rows[i]; 12236f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 12246f0a148fSBarry Smith } 12256f0a148fSBarry Smith if (diag) { 12266f0a148fSBarry Smith for (i=0; i<N; i++) { 12276f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1228260925b8SBarry Smith *slot = *diag; 12296f0a148fSBarry Smith } 12306f0a148fSBarry Smith } 1231260925b8SBarry Smith ISRestoreIndices(is,&rows); 12323a40ed3dSBarry Smith PetscFunctionReturn(0); 12336f0a148fSBarry Smith } 1234557bce09SLois Curfman McInnes 12354a2ae208SSatish Balay #undef __FUNCT__ 12364a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 123787828ca2SBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar **array) 123864e87e97SBarry Smith { 1239c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12403a40ed3dSBarry Smith 12413a40ed3dSBarry Smith PetscFunctionBegin; 124264e87e97SBarry Smith *array = mat->v; 12433a40ed3dSBarry Smith PetscFunctionReturn(0); 124464e87e97SBarry Smith } 12450754003eSLois Curfman McInnes 12464a2ae208SSatish Balay #undef __FUNCT__ 12474a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 124887828ca2SBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar **array) 1249ff14e315SSatish Balay { 12503a40ed3dSBarry Smith PetscFunctionBegin; 125109b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12523a40ed3dSBarry Smith PetscFunctionReturn(0); 1253ff14e315SSatish Balay } 12540754003eSLois Curfman McInnes 12554a2ae208SSatish Balay #undef __FUNCT__ 12564a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 12577b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12580754003eSLois Curfman McInnes { 1259c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1260273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 126187828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 12620754003eSLois Curfman McInnes Mat newmat; 12630754003eSLois Curfman McInnes 12643a40ed3dSBarry Smith PetscFunctionBegin; 126578b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 126678b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1267e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1268e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 12690754003eSLois Curfman McInnes 1270182d2002SSatish Balay /* Check submatrixcall */ 1271182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1272182d2002SSatish Balay int n_cols,n_rows; 1273182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 127429bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1275182d2002SSatish Balay newmat = *B; 1276182d2002SSatish Balay } else { 12770754003eSLois Curfman McInnes /* Create and fill new matrix */ 1278b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1279182d2002SSatish Balay } 1280182d2002SSatish Balay 1281182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1282182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1283182d2002SSatish Balay 1284182d2002SSatish Balay for (i=0; i<ncols; i++) { 1285182d2002SSatish Balay av = v + m*icol[i]; 1286182d2002SSatish Balay for (j=0; j<nrows; j++) { 1287182d2002SSatish Balay *bv++ = av[irow[j]]; 12880754003eSLois Curfman McInnes } 12890754003eSLois Curfman McInnes } 1290182d2002SSatish Balay 1291182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 12926d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12936d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12940754003eSLois Curfman McInnes 12950754003eSLois Curfman McInnes /* Free work space */ 129678b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 129778b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1298182d2002SSatish Balay *B = newmat; 12993a40ed3dSBarry Smith PetscFunctionReturn(0); 13000754003eSLois Curfman McInnes } 13010754003eSLois Curfman McInnes 13024a2ae208SSatish Balay #undef __FUNCT__ 13034a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 13047b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1305905e6a2fSBarry Smith { 1306905e6a2fSBarry Smith int ierr,i; 1307905e6a2fSBarry Smith 13083a40ed3dSBarry Smith PetscFunctionBegin; 1309905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1310b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1311905e6a2fSBarry Smith } 1312905e6a2fSBarry Smith 1313905e6a2fSBarry Smith for (i=0; i<n; i++) { 13146a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1315905e6a2fSBarry Smith } 13163a40ed3dSBarry Smith PetscFunctionReturn(0); 1317905e6a2fSBarry Smith } 1318905e6a2fSBarry Smith 13194a2ae208SSatish Balay #undef __FUNCT__ 13204a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1321cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 13224b0e389bSBarry Smith { 13234b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 13243a40ed3dSBarry Smith int ierr; 1325273d9f13SBarry Smith PetscTruth flag; 13263a40ed3dSBarry Smith 13273a40ed3dSBarry Smith PetscFunctionBegin; 1328273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr); 1329273d9f13SBarry Smith if (!flag) { 1330cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 13313a40ed3dSBarry Smith PetscFunctionReturn(0); 13323a40ed3dSBarry Smith } 1333273d9f13SBarry Smith if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 133487828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1335273d9f13SBarry Smith PetscFunctionReturn(0); 1336273d9f13SBarry Smith } 1337273d9f13SBarry Smith 13384a2ae208SSatish Balay #undef __FUNCT__ 13394a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1340273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1341273d9f13SBarry Smith { 1342273d9f13SBarry Smith int ierr; 1343273d9f13SBarry Smith 1344273d9f13SBarry Smith PetscFunctionBegin; 1345273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 13463a40ed3dSBarry Smith PetscFunctionReturn(0); 13474b0e389bSBarry Smith } 13484b0e389bSBarry Smith 1349289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1350a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1351905e6a2fSBarry Smith MatGetRow_SeqDense, 1352905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1353905e6a2fSBarry Smith MatMult_SeqDense, 1354905e6a2fSBarry Smith MatMultAdd_SeqDense, 13557c922b88SBarry Smith MatMultTranspose_SeqDense, 13567c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1357905e6a2fSBarry Smith MatSolve_SeqDense, 1358905e6a2fSBarry Smith MatSolveAdd_SeqDense, 13597c922b88SBarry Smith MatSolveTranspose_SeqDense, 13607c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1361905e6a2fSBarry Smith MatLUFactor_SeqDense, 1362905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1363ec8511deSBarry Smith MatRelax_SeqDense, 1364ec8511deSBarry Smith MatTranspose_SeqDense, 1365905e6a2fSBarry Smith MatGetInfo_SeqDense, 1366905e6a2fSBarry Smith MatEqual_SeqDense, 1367905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1368905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1369905e6a2fSBarry Smith MatNorm_SeqDense, 1370905e6a2fSBarry Smith 0, 1371905e6a2fSBarry Smith 0, 1372905e6a2fSBarry Smith 0, 1373905e6a2fSBarry Smith MatSetOption_SeqDense, 1374905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1375905e6a2fSBarry Smith MatZeroRows_SeqDense, 1376905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1377905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1378905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1379905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1380273d9f13SBarry Smith MatSetUpPreallocation_SeqDense, 1381273d9f13SBarry Smith 0, 1382905e6a2fSBarry Smith 0, 1383905e6a2fSBarry Smith MatGetArray_SeqDense, 1384905e6a2fSBarry Smith MatRestoreArray_SeqDense, 13855609ef8eSBarry Smith MatDuplicate_SeqDense, 1386a5ae1ecdSBarry Smith 0, 1387a5ae1ecdSBarry Smith 0, 1388a5ae1ecdSBarry Smith 0, 1389a5ae1ecdSBarry Smith 0, 1390a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1391a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1392a5ae1ecdSBarry Smith 0, 13934b0e389bSBarry Smith MatGetValues_SeqDense, 1394a5ae1ecdSBarry Smith MatCopy_SeqDense, 1395a5ae1ecdSBarry Smith 0, 1396a5ae1ecdSBarry Smith MatScale_SeqDense, 1397a5ae1ecdSBarry Smith 0, 1398a5ae1ecdSBarry Smith 0, 1399a5ae1ecdSBarry Smith 0, 1400a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1401a5ae1ecdSBarry Smith 0, 1402a5ae1ecdSBarry Smith 0, 1403a5ae1ecdSBarry Smith 0, 1404a5ae1ecdSBarry Smith 0, 1405a5ae1ecdSBarry Smith 0, 1406a5ae1ecdSBarry Smith 0, 1407a5ae1ecdSBarry Smith 0, 1408a5ae1ecdSBarry Smith 0, 1409a5ae1ecdSBarry Smith 0, 1410a5ae1ecdSBarry Smith 0, 1411e03a110bSBarry Smith MatDestroy_SeqDense, 1412e03a110bSBarry Smith MatView_SeqDense, 14138a124369SBarry Smith MatGetPetscMaps_Petsc}; 141490ace30eSBarry Smith 14154a2ae208SSatish Balay #undef __FUNCT__ 14164a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 14174b828684SBarry Smith /*@C 1418fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1419d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1420d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1421289bc588SBarry Smith 1422db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1423db81eaa0SLois Curfman McInnes 142420563c6bSBarry Smith Input Parameters: 1425db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 14260c775827SLois Curfman McInnes . m - number of rows 142718f449edSLois Curfman McInnes . n - number of columns 1428db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1429dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 143020563c6bSBarry Smith 143120563c6bSBarry Smith Output Parameter: 143244cd7ae7SLois Curfman McInnes . A - the matrix 143320563c6bSBarry Smith 1434b259b22eSLois Curfman McInnes Notes: 143518f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 143618f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1437b4fd4287SBarry Smith set data=PETSC_NULL. 143818f449edSLois Curfman McInnes 1439027ccd11SLois Curfman McInnes Level: intermediate 1440027ccd11SLois Curfman McInnes 1441dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1442d65003e9SLois Curfman McInnes 1443db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 144420563c6bSBarry Smith @*/ 144587828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1446289bc588SBarry Smith { 1447273d9f13SBarry Smith int ierr; 14483b2fbd54SBarry Smith 14493a40ed3dSBarry Smith PetscFunctionBegin; 1450273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1451273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1452273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1453273d9f13SBarry Smith PetscFunctionReturn(0); 1454273d9f13SBarry Smith } 1455273d9f13SBarry Smith 14564a2ae208SSatish Balay #undef __FUNCT__ 14574a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1458273d9f13SBarry Smith /*@C 1459273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1460273d9f13SBarry Smith 1461273d9f13SBarry Smith Collective on MPI_Comm 1462273d9f13SBarry Smith 1463273d9f13SBarry Smith Input Parameters: 1464273d9f13SBarry Smith + A - the matrix 1465273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1466273d9f13SBarry Smith 1467273d9f13SBarry Smith Notes: 1468273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1469273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1470273d9f13SBarry Smith set data=PETSC_NULL. 1471273d9f13SBarry Smith 1472273d9f13SBarry Smith Level: intermediate 1473273d9f13SBarry Smith 1474273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1475273d9f13SBarry Smith 1476273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1477273d9f13SBarry Smith @*/ 147887828ca2SBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data) 1479273d9f13SBarry Smith { 1480273d9f13SBarry Smith Mat_SeqDense *b; 1481273d9f13SBarry Smith int ierr; 1482273d9f13SBarry Smith PetscTruth flg2; 1483273d9f13SBarry Smith 1484273d9f13SBarry Smith PetscFunctionBegin; 1485273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr); 1486273d9f13SBarry Smith if (!flg2) PetscFunctionReturn(0); 1487273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1488273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1489273d9f13SBarry Smith if (!data) { 149087828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 149187828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1492273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 149387828ca2SBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1494273d9f13SBarry Smith } else { /* user-allocated storage */ 1495273d9f13SBarry Smith b->v = data; 1496273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1497273d9f13SBarry Smith } 1498273d9f13SBarry Smith PetscFunctionReturn(0); 1499273d9f13SBarry Smith } 1500273d9f13SBarry Smith 1501273d9f13SBarry Smith EXTERN_C_BEGIN 15024a2ae208SSatish Balay #undef __FUNCT__ 15034a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1504273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1505273d9f13SBarry Smith { 1506273d9f13SBarry Smith Mat_SeqDense *b; 1507273d9f13SBarry Smith int ierr,size; 1508273d9f13SBarry Smith 1509273d9f13SBarry Smith PetscFunctionBegin; 1510273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 151129bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 151255659b69SBarry Smith 1513273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1514273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1515273d9f13SBarry Smith 1516b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1517549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1518549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 151944cd7ae7SLois Curfman McInnes B->factor = 0; 152090f02eecSBarry Smith B->mapping = 0; 1521b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 152244cd7ae7SLois Curfman McInnes B->data = (void*)b; 152318f449edSLois Curfman McInnes 15248a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 15258a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1526a5ae1ecdSBarry Smith 152744cd7ae7SLois Curfman McInnes b->pivots = 0; 1528273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1529273d9f13SBarry Smith b->v = 0; 15304e220ebcSLois Curfman McInnes 15313a40ed3dSBarry Smith PetscFunctionReturn(0); 1532289bc588SBarry Smith } 1533273d9f13SBarry Smith EXTERN_C_END 15343b2fbd54SBarry Smith 15353b2fbd54SBarry Smith 15363b2fbd54SBarry Smith 15373b2fbd54SBarry Smith 1538