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; 17*1b807ce4Svictorle if (x->lda>X->m || y->lda>Y->m) SETERRQ(1,"Can not handle LDA"); 181987afe7SBarry Smith BLaxpy_(&N,alpha,x->v,&one,y->v,&one); 19b0a32e0cSBarry Smith PetscLogFlops(2*N-1); 203a40ed3dSBarry Smith PetscFunctionReturn(0); 211987afe7SBarry Smith } 221987afe7SBarry Smith 234a2ae208SSatish Balay #undef __FUNCT__ 244a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 258f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 26289bc588SBarry Smith { 274e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 28273d9f13SBarry Smith int i,N = A->m*A->n,count = 0; 2987828ca2SBarry Smith PetscScalar *v = a->v; 303a40ed3dSBarry Smith 313a40ed3dSBarry Smith PetscFunctionBegin; 32289bc588SBarry Smith for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 334e220ebcSLois Curfman McInnes 34273d9f13SBarry Smith info->rows_global = (double)A->m; 35273d9f13SBarry Smith info->columns_global = (double)A->n; 36273d9f13SBarry Smith info->rows_local = (double)A->m; 37273d9f13SBarry Smith info->columns_local = (double)A->n; 384e220ebcSLois Curfman McInnes info->block_size = 1.0; 394e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 404e220ebcSLois Curfman McInnes info->nz_used = (double)count; 414e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 424e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 434e220ebcSLois Curfman McInnes info->mallocs = 0; 444e220ebcSLois Curfman McInnes info->memory = A->mem; 454e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 464e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 474e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 484e220ebcSLois Curfman McInnes 493a40ed3dSBarry Smith PetscFunctionReturn(0); 50289bc588SBarry Smith } 51289bc588SBarry Smith 524a2ae208SSatish Balay #undef __FUNCT__ 534a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 5487828ca2SBarry Smith int MatScale_SeqDense(PetscScalar *alpha,Mat A) 5580cd9d93SLois Curfman McInnes { 56273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 5780cd9d93SLois Curfman McInnes int one = 1,nz; 5880cd9d93SLois Curfman McInnes 593a40ed3dSBarry Smith PetscFunctionBegin; 60*1b807ce4Svictorle if (a->lda>A->m) SETERRQ(1,"Can not handle LDA"); 61273d9f13SBarry Smith nz = A->m*A->n; 6280cd9d93SLois Curfman McInnes BLscal_(&nz,alpha,a->v,&one); 63b0a32e0cSBarry Smith PetscLogFlops(nz); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 6580cd9d93SLois Curfman McInnes } 6680cd9d93SLois Curfman McInnes 67289bc588SBarry Smith /* ---------------------------------------------------------------*/ 686831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 696831982aSBarry Smith rather than put it in the Mat->row slot.*/ 704a2ae208SSatish Balay #undef __FUNCT__ 714a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 72b380c88cSHong Zhang int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 73289bc588SBarry Smith { 74c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 75b0a32e0cSBarry Smith int info,ierr; 7655659b69SBarry Smith 773a40ed3dSBarry Smith PetscFunctionBegin; 78289bc588SBarry Smith if (!mat->pivots) { 79b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr); 80b0a32e0cSBarry Smith PetscLogObjectMemory(A,A->m*sizeof(int)); 81289bc588SBarry Smith } 82f1af5d2fSBarry Smith A->factor = FACTOR_LU; 83273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 84ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRF) 85ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable."); 86ae7cfcebSSatish Balay #else 87*1b807ce4Svictorle LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info); 8829bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 8929bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 902ffef20aSBarry Smith #endif 91b0a32e0cSBarry Smith PetscLogFlops((2*A->n*A->n*A->n)/3); 923a40ed3dSBarry Smith PetscFunctionReturn(0); 93289bc588SBarry Smith } 946ee01492SSatish Balay 954a2ae208SSatish Balay #undef __FUNCT__ 964a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 975609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9802cad45dSBarry Smith { 9902cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 10002cad45dSBarry Smith int ierr; 10102cad45dSBarry Smith Mat newi; 10202cad45dSBarry Smith 1033a40ed3dSBarry Smith PetscFunctionBegin; 104*1b807ce4Svictorle if (mat->lda>A->m) SETERRQ(1,"Can not handle LDA"); 105273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr); 10602cad45dSBarry Smith l = (Mat_SeqDense*)newi->data; 1075609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 10887828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 10902cad45dSBarry Smith } 11002cad45dSBarry Smith newi->assembled = PETSC_TRUE; 11102cad45dSBarry Smith *newmat = newi; 1123a40ed3dSBarry Smith PetscFunctionReturn(0); 11302cad45dSBarry Smith } 11402cad45dSBarry Smith 1154a2ae208SSatish Balay #undef __FUNCT__ 1164a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 117b380c88cSHong Zhang int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 118289bc588SBarry Smith { 1193a40ed3dSBarry Smith int ierr; 1203a40ed3dSBarry Smith 1213a40ed3dSBarry Smith PetscFunctionBegin; 1225609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1233a40ed3dSBarry Smith PetscFunctionReturn(0); 124289bc588SBarry Smith } 1256ee01492SSatish Balay 1264a2ae208SSatish Balay #undef __FUNCT__ 1274a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1288f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 129289bc588SBarry Smith { 13002cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 131*1b807ce4Svictorle int lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr; 1323a40ed3dSBarry Smith 1333a40ed3dSBarry Smith PetscFunctionBegin; 13402cad45dSBarry Smith /* copy the numerical values */ 135*1b807ce4Svictorle if (lda1>m || lda2>m ) { 136*1b807ce4Svictorle for (j=0; j<n; j++) { 137*1b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar)); CHKERRQ(ierr); 138*1b807ce4Svictorle } 139*1b807ce4Svictorle } else { 14087828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 141*1b807ce4Svictorle } 14202cad45dSBarry Smith (*fact)->factor = 0; 143e03a110bSBarry Smith ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr); 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 145289bc588SBarry Smith } 1466ee01492SSatish Balay 1474a2ae208SSatish Balay #undef __FUNCT__ 1484a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 14915e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 150289bc588SBarry Smith { 1513a40ed3dSBarry Smith int ierr; 1523a40ed3dSBarry Smith 1533a40ed3dSBarry Smith PetscFunctionBegin; 1543a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156289bc588SBarry Smith } 1576ee01492SSatish Balay 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 16015e8a5b3SHong Zhang int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 161289bc588SBarry Smith { 162c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 163606d414cSSatish Balay int info,ierr; 16455659b69SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 166752f5567SLois Curfman McInnes if (mat->pivots) { 167606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 168b0a32e0cSBarry Smith PetscLogObjectMemory(A,-A->m*sizeof(int)); 169752f5567SLois Curfman McInnes mat->pivots = 0; 170752f5567SLois Curfman McInnes } 171273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 172ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRF) 173ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable."); 174ae7cfcebSSatish Balay #else 175*1b807ce4Svictorle LApotrf_("L",&A->n,mat->v,&mat->lda,&info); 17629bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 1772ffef20aSBarry Smith #endif 178c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 179b0a32e0cSBarry Smith PetscLogFlops((A->n*A->n*A->n)/3); 1803a40ed3dSBarry Smith PetscFunctionReturn(0); 181289bc588SBarry Smith } 182289bc588SBarry Smith 1834a2ae208SSatish Balay #undef __FUNCT__ 1840b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 1850b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 1860b4b3355SBarry Smith { 1870b4b3355SBarry Smith int ierr; 18815e8a5b3SHong Zhang MatFactorInfo info; 1890b4b3355SBarry Smith 1900b4b3355SBarry Smith PetscFunctionBegin; 19115e8a5b3SHong Zhang info.fill = 1.0; 19215e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 1930b4b3355SBarry Smith PetscFunctionReturn(0); 1940b4b3355SBarry Smith } 1950b4b3355SBarry Smith 1960b4b3355SBarry Smith #undef __FUNCT__ 1974a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 1988f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 199289bc588SBarry Smith { 200c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 201a57ad546SLois Curfman McInnes int one = 1,info,ierr; 20287828ca2SBarry Smith PetscScalar *x,*y; 20367e560aaSBarry Smith 2043a40ed3dSBarry Smith PetscFunctionBegin; 205273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2067a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2077a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 20887828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 209c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 210ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 211ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 212ae7cfcebSSatish Balay #else 213*1b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2142ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 215ae7cfcebSSatish Balay #endif 2167a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 217ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 218ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 219ae7cfcebSSatish Balay #else 220*1b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2212ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 222ae7cfcebSSatish Balay #endif 223289bc588SBarry Smith } 22429bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2257a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2267a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 227b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229289bc588SBarry Smith } 2306ee01492SSatish Balay 2314a2ae208SSatish Balay #undef __FUNCT__ 2324a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 2337c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 234da3a660dSBarry Smith { 235c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2367a97a34bSBarry Smith int ierr,one = 1,info; 23787828ca2SBarry Smith PetscScalar *x,*y; 23867e560aaSBarry Smith 2393a40ed3dSBarry Smith PetscFunctionBegin; 240273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2417a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2427a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 24387828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 244752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 245da3a660dSBarry Smith if (mat->pivots) { 246ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 247ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 248ae7cfcebSSatish Balay #else 249*1b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2502ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 251ae7cfcebSSatish Balay #endif 2527a97a34bSBarry Smith } else { 253ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 254ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 255ae7cfcebSSatish Balay #else 256*1b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2572ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 258ae7cfcebSSatish Balay #endif 259da3a660dSBarry Smith } 2607a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2617a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 262b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 264da3a660dSBarry Smith } 2656ee01492SSatish Balay 2664a2ae208SSatish Balay #undef __FUNCT__ 2674a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 2688f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 269da3a660dSBarry Smith { 270c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2717a97a34bSBarry Smith int ierr,one = 1,info; 27287828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 273da3a660dSBarry Smith Vec tmp = 0; 27467e560aaSBarry Smith 2753a40ed3dSBarry Smith PetscFunctionBegin; 2767a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2777a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 278273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 279da3a660dSBarry Smith if (yy == zz) { 28078b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 281b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 28278b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 283da3a660dSBarry Smith } 28487828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 285752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 286da3a660dSBarry Smith if (mat->pivots) { 287ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 288ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 289ae7cfcebSSatish Balay #else 290*1b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2912ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 292ae7cfcebSSatish Balay #endif 293a8c6a408SBarry Smith } else { 294ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 295ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 296ae7cfcebSSatish Balay #else 297*1b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2982ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 299ae7cfcebSSatish Balay #endif 300da3a660dSBarry Smith } 301a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 302a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 3037a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3047a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 305b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 307da3a660dSBarry Smith } 30867e560aaSBarry Smith 3094a2ae208SSatish Balay #undef __FUNCT__ 3104a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 3117c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 312da3a660dSBarry Smith { 313c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3146abc6512SBarry Smith int one = 1,info,ierr; 31587828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 316da3a660dSBarry Smith Vec tmp; 31767e560aaSBarry Smith 3183a40ed3dSBarry Smith PetscFunctionBegin; 319273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3207a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3217a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 322da3a660dSBarry Smith if (yy == zz) { 32378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 324b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 32578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 326da3a660dSBarry Smith } 32787828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 328752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 329da3a660dSBarry Smith if (mat->pivots) { 330ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 331ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 332ae7cfcebSSatish Balay #else 333*1b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 3342ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 335ae7cfcebSSatish Balay #endif 3363a40ed3dSBarry Smith } else { 337ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 338ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 339ae7cfcebSSatish Balay #else 340*1b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 3412ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 342ae7cfcebSSatish Balay #endif 343da3a660dSBarry Smith } 34490f02eecSBarry Smith if (tmp) { 34590f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 34690f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3473a40ed3dSBarry Smith } else { 34890f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 34990f02eecSBarry Smith } 3507a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3517a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 352b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3533a40ed3dSBarry Smith PetscFunctionReturn(0); 354da3a660dSBarry Smith } 355289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3564a2ae208SSatish Balay #undef __FUNCT__ 3574a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 358329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 359c14dc6b6SHong Zhang PetscReal shift,int its,int lits,Vec xx) 360289bc588SBarry Smith { 361c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 36287828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 363273d9f13SBarry Smith int ierr,m = A->m,i; 364aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 365bc1b551cSSatish Balay int o = 1; 366bc1b551cSSatish Balay #endif 367289bc588SBarry Smith 3683a40ed3dSBarry Smith PetscFunctionBegin; 369289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3703a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 371bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 372289bc588SBarry Smith } 3737a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3747a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 375b965ef7fSBarry Smith its = its*lits; 37691723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 377289bc588SBarry Smith while (its--) { 378289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 379289bc588SBarry Smith for (i=0; i<m; 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; 38487828ca2SBarry Smith PetscScalar 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 if (flag & SOR_BACKWARD_SWEEP) { 396289bc588SBarry Smith for (i=m-1; i>=0; i--) { 397aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 398f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 399f1747703SBarry Smith not happy about returning a double complex */ 400f1747703SBarry Smith int _i; 40187828ca2SBarry Smith PetscScalar sum = b[i]; 402f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4033f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 404f1747703SBarry Smith } 405f1747703SBarry Smith xt = sum; 406f1747703SBarry Smith #else 407289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 408f1747703SBarry Smith #endif 40955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 410289bc588SBarry Smith } 411289bc588SBarry Smith } 412289bc588SBarry Smith } 413600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4147a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4153a40ed3dSBarry Smith PetscFunctionReturn(0); 416289bc588SBarry Smith } 417289bc588SBarry Smith 418289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 4217c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 422289bc588SBarry Smith { 423c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 42487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 425ea709b57SSatish Balay int ierr,_One=1; 426ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4273a40ed3dSBarry Smith 4283a40ed3dSBarry Smith PetscFunctionBegin; 429273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4307a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4317a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 432*1b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4337a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4347a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 435b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->n); 4363a40ed3dSBarry Smith PetscFunctionReturn(0); 437289bc588SBarry Smith } 4386ee01492SSatish Balay 4394a2ae208SSatish Balay #undef __FUNCT__ 4404a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 44144cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 442289bc588SBarry Smith { 443c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 44487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 445329f5518SBarry Smith int ierr,_One=1; 4463a40ed3dSBarry Smith 4473a40ed3dSBarry Smith PetscFunctionBegin; 448273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4497a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4507a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 451*1b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4527a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4537a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 454b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->m); 4553a40ed3dSBarry Smith PetscFunctionReturn(0); 456289bc588SBarry Smith } 4576ee01492SSatish Balay 4584a2ae208SSatish Balay #undef __FUNCT__ 4594a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 46044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 461289bc588SBarry Smith { 462c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 46387828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 464329f5518SBarry Smith int ierr,_One=1; 4653a40ed3dSBarry Smith 4663a40ed3dSBarry Smith PetscFunctionBegin; 467273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 468600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4697a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4707a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 471*1b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 4727a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4737a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 474b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4753a40ed3dSBarry Smith PetscFunctionReturn(0); 476289bc588SBarry Smith } 4776ee01492SSatish Balay 4784a2ae208SSatish Balay #undef __FUNCT__ 4794a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 4807c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 481289bc588SBarry Smith { 482c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 48387828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 4847a97a34bSBarry Smith int ierr,_One=1; 48587828ca2SBarry Smith PetscScalar _DOne=1.0; 4863a40ed3dSBarry Smith 4873a40ed3dSBarry Smith PetscFunctionBegin; 488273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 489600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4907a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4917a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 492*1b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 4937a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4947a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 495b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4963a40ed3dSBarry Smith PetscFunctionReturn(0); 497289bc588SBarry Smith } 498289bc588SBarry Smith 499289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5004a2ae208SSatish Balay #undef __FUNCT__ 5014a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 50287828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 503289bc588SBarry Smith { 504c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 50587828ca2SBarry Smith PetscScalar *v; 506b0a32e0cSBarry Smith int i,ierr; 50767e560aaSBarry Smith 5083a40ed3dSBarry Smith PetscFunctionBegin; 509273d9f13SBarry Smith *ncols = A->n; 510289bc588SBarry Smith if (cols) { 511b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 512273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 513289bc588SBarry Smith } 514289bc588SBarry Smith if (vals) { 51587828ca2SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 516289bc588SBarry Smith v = mat->v + row; 517*1b807ce4Svictorle for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;} 518289bc588SBarry Smith } 5193a40ed3dSBarry Smith PetscFunctionReturn(0); 520289bc588SBarry Smith } 5216ee01492SSatish Balay 5224a2ae208SSatish Balay #undef __FUNCT__ 5234a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 52487828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 525289bc588SBarry Smith { 526606d414cSSatish Balay int ierr; 527606d414cSSatish Balay PetscFunctionBegin; 528606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 529606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5303a40ed3dSBarry Smith PetscFunctionReturn(0); 531289bc588SBarry Smith } 532289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5334a2ae208SSatish Balay #undef __FUNCT__ 5344a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 5358f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 53687828ca2SBarry Smith int *indexn,PetscScalar *v,InsertMode addv) 537289bc588SBarry Smith { 538c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 539289bc588SBarry Smith int i,j; 540d6dfbf8fSBarry Smith 5413a40ed3dSBarry Smith PetscFunctionBegin; 542289bc588SBarry Smith if (!mat->roworiented) { 543dbb450caSBarry Smith if (addv == INSERT_VALUES) { 544289bc588SBarry Smith for (j=0; j<n; j++) { 5455ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 546aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 54729bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 54858804f6dSBarry Smith #endif 549289bc588SBarry Smith for (i=0; i<m; i++) { 5505ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 551aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 55229bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 55358804f6dSBarry Smith #endif 554*1b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] = *v++; 555289bc588SBarry Smith } 556289bc588SBarry Smith } 5573a40ed3dSBarry Smith } else { 558289bc588SBarry Smith for (j=0; j<n; j++) { 5595ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 560aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 56129bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56258804f6dSBarry Smith #endif 563289bc588SBarry Smith for (i=0; i<m; i++) { 5645ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 565aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 56629bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56758804f6dSBarry Smith #endif 568*1b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] += *v++; 569289bc588SBarry Smith } 570289bc588SBarry Smith } 571289bc588SBarry Smith } 5723a40ed3dSBarry Smith } else { 573dbb450caSBarry Smith if (addv == INSERT_VALUES) { 574e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 5755ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 576aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57729bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 57858804f6dSBarry Smith #endif 579e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 5805ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 581aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 58229bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 58358804f6dSBarry Smith #endif 584*1b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] = *v++; 585e8d4e0b9SBarry Smith } 586e8d4e0b9SBarry Smith } 5873a40ed3dSBarry Smith } else { 588289bc588SBarry Smith for (i=0; i<m; i++) { 5895ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 590aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 59129bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 59258804f6dSBarry Smith #endif 593289bc588SBarry Smith for (j=0; j<n; j++) { 5945ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 595aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 59629bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 59758804f6dSBarry Smith #endif 598*1b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] += *v++; 599289bc588SBarry Smith } 600289bc588SBarry Smith } 601289bc588SBarry Smith } 602e8d4e0b9SBarry Smith } 6033a40ed3dSBarry Smith PetscFunctionReturn(0); 604289bc588SBarry Smith } 605e8d4e0b9SBarry Smith 6064a2ae208SSatish Balay #undef __FUNCT__ 6074a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 60887828ca2SBarry Smith int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v) 609ae80bb75SLois Curfman McInnes { 610ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 611ae80bb75SLois Curfman McInnes int i,j; 61287828ca2SBarry Smith PetscScalar *vpt = v; 613ae80bb75SLois Curfman McInnes 6143a40ed3dSBarry Smith PetscFunctionBegin; 615ae80bb75SLois Curfman McInnes /* row-oriented output */ 616ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 617ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 618*1b807ce4Svictorle *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 619ae80bb75SLois Curfman McInnes } 620ae80bb75SLois Curfman McInnes } 6213a40ed3dSBarry Smith PetscFunctionReturn(0); 622ae80bb75SLois Curfman McInnes } 623ae80bb75SLois Curfman McInnes 624289bc588SBarry Smith /* -----------------------------------------------------------------*/ 625289bc588SBarry Smith 626e090d566SSatish Balay #include "petscsys.h" 62788e32edaSLois Curfman McInnes 628273d9f13SBarry Smith EXTERN_C_BEGIN 6294a2ae208SSatish Balay #undef __FUNCT__ 6304a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 631b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A) 63288e32edaSLois Curfman McInnes { 63388e32edaSLois Curfman McInnes Mat_SeqDense *a; 63488e32edaSLois Curfman McInnes Mat B; 63588e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 63688e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 63787828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 63819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 63988e32edaSLois Curfman McInnes 6403a40ed3dSBarry Smith PetscFunctionBegin; 641d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 64229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 643b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6440752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 645552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 64688e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 64788e32edaSLois Curfman McInnes 64890ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 64990ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 65090ace30eSBarry Smith B = *A; 65190ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6524a41a4d3SSatish Balay v = a->v; 6534a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6544a41a4d3SSatish Balay from row major to column major */ 65587828ca2SBarry Smith ierr = PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 65690ace30eSBarry Smith /* read in nonzero values */ 6574a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6584a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6594a41a4d3SSatish Balay for (j=0; j<N; j++) { 6604a41a4d3SSatish Balay for (i=0; i<M; i++) { 6614a41a4d3SSatish Balay *v++ =w[i*N+j]; 6624a41a4d3SSatish Balay } 6634a41a4d3SSatish Balay } 664606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6656d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6666d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 66790ace30eSBarry Smith } else { 66888e32edaSLois Curfman McInnes /* read row lengths */ 669b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 6700752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 67188e32edaSLois Curfman McInnes 67288e32edaSLois Curfman McInnes /* create our matrix */ 673b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 67488e32edaSLois Curfman McInnes B = *A; 67588e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 67688e32edaSLois Curfman McInnes v = a->v; 67788e32edaSLois Curfman McInnes 67888e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 679b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 680b0a32e0cSBarry Smith cols = scols; 6810752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 68287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 683b0a32e0cSBarry Smith vals = svals; 6840752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 68588e32edaSLois Curfman McInnes 68688e32edaSLois Curfman McInnes /* insert into matrix */ 68788e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 68888e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 68988e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 69088e32edaSLois Curfman McInnes } 691606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 692606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 693606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 69488e32edaSLois Curfman McInnes 6956d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6966d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69790ace30eSBarry Smith } 6983a40ed3dSBarry Smith PetscFunctionReturn(0); 69988e32edaSLois Curfman McInnes } 700273d9f13SBarry Smith EXTERN_C_END 70188e32edaSLois Curfman McInnes 702e090d566SSatish Balay #include "petscsys.h" 703932b0c3eSLois Curfman McInnes 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 706b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 707289bc588SBarry Smith { 708932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 709fb9695e5SSatish Balay int ierr,i,j; 710fb9695e5SSatish Balay char *name; 71187828ca2SBarry Smith PetscScalar *v; 712f3ef73ceSBarry Smith PetscViewerFormat format; 713932b0c3eSLois Curfman McInnes 7143a40ed3dSBarry Smith PetscFunctionBegin; 715435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 716b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 717456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7183a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 719fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 720b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 721273d9f13SBarry Smith for (i=0; i<A->m; i++) { 72244cd7ae7SLois Curfman McInnes v = a->v + i; 723b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 724273d9f13SBarry Smith for (j=0; j<A->n; j++) { 725aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 726329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 72762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 728329f5518SBarry Smith } else if (PetscRealPart(*v)) { 72962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7306831982aSBarry Smith } 73180cd9d93SLois Curfman McInnes #else 7326831982aSBarry Smith if (*v) { 73362b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 7346831982aSBarry Smith } 73580cd9d93SLois Curfman McInnes #endif 736*1b807ce4Svictorle v += a->lda; 73780cd9d93SLois Curfman McInnes } 738b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 73980cd9d93SLois Curfman McInnes } 740b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7413a40ed3dSBarry Smith } else { 742b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 743aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 744ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 74547989497SBarry Smith /* determine if matrix has all real values */ 74647989497SBarry Smith v = a->v; 747*1b807ce4Svictorle for (i=0; i<A->m; i++) { 748*1b807ce4Svictorle v = a->v + i; 749*1b807ce4Svictorle for (j=0; j<A->n; j++) { 750ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 751*1b807ce4Svictorle v += a->lda; 752*1b807ce4Svictorle } 75347989497SBarry Smith } 75447989497SBarry Smith #endif 755fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7563a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 757b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 758fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 759fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 760ffac6cdbSBarry Smith } 761ffac6cdbSBarry Smith 762273d9f13SBarry Smith for (i=0; i<A->m; i++) { 763932b0c3eSLois Curfman McInnes v = a->v + i; 764273d9f13SBarry Smith for (j=0; j<A->n; j++) { 765aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 76647989497SBarry Smith if (allreal) { 767*1b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 76847989497SBarry Smith } else { 769*1b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 77047989497SBarry Smith } 771289bc588SBarry Smith #else 772*1b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 773289bc588SBarry Smith #endif 774*1b807ce4Svictorle v += a->lda; 775289bc588SBarry Smith } 776b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 777289bc588SBarry Smith } 778fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 779b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 780ffac6cdbSBarry Smith } 781b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 782da3a660dSBarry Smith } 783b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7843a40ed3dSBarry Smith PetscFunctionReturn(0); 785289bc588SBarry Smith } 786289bc588SBarry Smith 7874a2ae208SSatish Balay #undef __FUNCT__ 7884a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 789b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 790932b0c3eSLois Curfman McInnes { 791932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 792273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 79387828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 794f3ef73ceSBarry Smith PetscViewerFormat format; 795932b0c3eSLois Curfman McInnes 7963a40ed3dSBarry Smith PetscFunctionBegin; 797b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 79890ace30eSBarry Smith 799b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 800fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 80190ace30eSBarry Smith /* store the matrix as a dense matrix */ 802b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 803552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 80490ace30eSBarry Smith col_lens[1] = m; 80590ace30eSBarry Smith col_lens[2] = n; 80690ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8070752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 808606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 80990ace30eSBarry Smith 81090ace30eSBarry Smith /* write out matrix, by rows */ 81187828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 81290ace30eSBarry Smith v = a->v; 81390ace30eSBarry Smith for (i=0; i<m; i++) { 81490ace30eSBarry Smith for (j=0; j<n; j++) { 81590ace30eSBarry Smith vals[i + j*m] = *v++; 81690ace30eSBarry Smith } 81790ace30eSBarry Smith } 8180752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 819606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 82090ace30eSBarry Smith } else { 821b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 822552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 823932b0c3eSLois Curfman McInnes col_lens[1] = m; 824932b0c3eSLois Curfman McInnes col_lens[2] = n; 825932b0c3eSLois Curfman McInnes col_lens[3] = nz; 826932b0c3eSLois Curfman McInnes 827932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 828932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8290752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 830932b0c3eSLois Curfman McInnes 831932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 832932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 833932b0c3eSLois Curfman McInnes ict = 0; 834932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 835932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 836932b0c3eSLois Curfman McInnes } 8370752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 838606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 839932b0c3eSLois Curfman McInnes 840932b0c3eSLois Curfman McInnes /* store nonzero values */ 84187828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 842932b0c3eSLois Curfman McInnes ict = 0; 843932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 844932b0c3eSLois Curfman McInnes v = a->v + i; 845932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 846*1b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 847932b0c3eSLois Curfman McInnes } 848932b0c3eSLois Curfman McInnes } 8490752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 850606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 85190ace30eSBarry Smith } 8523a40ed3dSBarry Smith PetscFunctionReturn(0); 853932b0c3eSLois Curfman McInnes } 854932b0c3eSLois Curfman McInnes 8554a2ae208SSatish Balay #undef __FUNCT__ 8564a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 857b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 858f1af5d2fSBarry Smith { 859f1af5d2fSBarry Smith Mat A = (Mat) Aa; 860f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 861fb9695e5SSatish Balay int m = A->m,n = A->n,color,i,j,ierr; 86287828ca2SBarry Smith PetscScalar *v = a->v; 863b0a32e0cSBarry Smith PetscViewer viewer; 864b0a32e0cSBarry Smith PetscDraw popup; 865329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 866f3ef73ceSBarry Smith PetscViewerFormat format; 867f1af5d2fSBarry Smith 868f1af5d2fSBarry Smith PetscFunctionBegin; 869f1af5d2fSBarry Smith 870f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 871b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 872b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 873f1af5d2fSBarry Smith 874f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 875fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 876f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 877b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 878f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 879f1af5d2fSBarry Smith x_l = j; 880f1af5d2fSBarry Smith x_r = x_l + 1.0; 881f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 882f1af5d2fSBarry Smith y_l = m - i - 1.0; 883f1af5d2fSBarry Smith y_r = y_l + 1.0; 884f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 885329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 886b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 887329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 888b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 889f1af5d2fSBarry Smith } else { 890f1af5d2fSBarry Smith continue; 891f1af5d2fSBarry Smith } 892f1af5d2fSBarry Smith #else 893f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 894b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 895f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 896b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 897f1af5d2fSBarry Smith } else { 898f1af5d2fSBarry Smith continue; 899f1af5d2fSBarry Smith } 900f1af5d2fSBarry Smith #endif 901b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 902f1af5d2fSBarry Smith } 903f1af5d2fSBarry Smith } 904f1af5d2fSBarry Smith } else { 905f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 906f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 907f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 908f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 909f1af5d2fSBarry Smith } 910b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 911b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 912b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 913f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 914f1af5d2fSBarry Smith x_l = j; 915f1af5d2fSBarry Smith x_r = x_l + 1.0; 916f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 917f1af5d2fSBarry Smith y_l = m - i - 1.0; 918f1af5d2fSBarry Smith y_r = y_l + 1.0; 919b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 920b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 921f1af5d2fSBarry Smith } 922f1af5d2fSBarry Smith } 923f1af5d2fSBarry Smith } 924f1af5d2fSBarry Smith PetscFunctionReturn(0); 925f1af5d2fSBarry Smith } 926f1af5d2fSBarry Smith 9274a2ae208SSatish Balay #undef __FUNCT__ 9284a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 929b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 930f1af5d2fSBarry Smith { 931b0a32e0cSBarry Smith PetscDraw draw; 932f1af5d2fSBarry Smith PetscTruth isnull; 933329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 934f1af5d2fSBarry Smith int ierr; 935f1af5d2fSBarry Smith 936f1af5d2fSBarry Smith PetscFunctionBegin; 937b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 938b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 939f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 940f1af5d2fSBarry Smith 941f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 942273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 943f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 944b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 945b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 946f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 947f1af5d2fSBarry Smith PetscFunctionReturn(0); 948f1af5d2fSBarry Smith } 949f1af5d2fSBarry Smith 9504a2ae208SSatish Balay #undef __FUNCT__ 9514a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 952b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer) 953932b0c3eSLois Curfman McInnes { 954932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 955bcd2baecSBarry Smith int ierr; 956f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 957932b0c3eSLois Curfman McInnes 9583a40ed3dSBarry Smith PetscFunctionBegin; 959b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 960b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 961fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 962fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 9630f5bd95cSBarry Smith 9640f5bd95cSBarry Smith if (issocket) { 965*1b807ce4Svictorle if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA"); 966b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9670f5bd95cSBarry Smith } else if (isascii) { 9683a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9690f5bd95cSBarry Smith } else if (isbinary) { 9703a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 971f1af5d2fSBarry Smith } else if (isdraw) { 972f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9735cd90555SBarry Smith } else { 97429bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 975932b0c3eSLois Curfman McInnes } 9763a40ed3dSBarry Smith PetscFunctionReturn(0); 977932b0c3eSLois Curfman McInnes } 978289bc588SBarry Smith 9794a2ae208SSatish Balay #undef __FUNCT__ 9804a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 981c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 982289bc588SBarry Smith { 983ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 98490f02eecSBarry Smith int ierr; 98590f02eecSBarry Smith 9863a40ed3dSBarry Smith PetscFunctionBegin; 987aa482453SBarry Smith #if defined(PETSC_USE_LOG) 988b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 989a5a9c739SBarry Smith #endif 990606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 991606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 992606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 9933a40ed3dSBarry Smith PetscFunctionReturn(0); 994289bc588SBarry Smith } 995289bc588SBarry Smith 9964a2ae208SSatish Balay #undef __FUNCT__ 9974a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 9988f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 999289bc588SBarry Smith { 1000c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1001*1b807ce4Svictorle int k,j,m,n,M,ierr; 100287828ca2SBarry Smith PetscScalar *v,tmp; 100348b35521SBarry Smith 10043a40ed3dSBarry Smith PetscFunctionBegin; 1005*1b807ce4Svictorle v = mat->v; m = A->m; M = mat->lda; n = A->n; 10067c922b88SBarry Smith if (!matout) { /* in place transpose */ 1007d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 100887828ca2SBarry Smith PetscScalar *w; 1009*1b807ce4Svictorle if (mat->lda>A->m) SETERRQ(1,"Can not handle LDA"); 101087828ca2SBarry Smith ierr = PetscMalloc((m+1)*(n+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 1011d64ed03dSBarry Smith for (j=0; j<m; j++) { 1012d64ed03dSBarry Smith for (k=0; k<n; k++) { 1013*1b807ce4Svictorle w[k + j*n] = v[j + k*M]; 1014d64ed03dSBarry Smith } 1015d64ed03dSBarry Smith } 101687828ca2SBarry Smith ierr = PetscMemcpy(v,w,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1017606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 1018d64ed03dSBarry Smith } else { 1019d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1020289bc588SBarry Smith for (k=0; k<j; k++) { 1021*1b807ce4Svictorle tmp = v[j + k*M]; 1022*1b807ce4Svictorle v[j + k*M] = v[k + j*M]; 1023*1b807ce4Svictorle v[k + j*M] = tmp; 1024289bc588SBarry Smith } 1025289bc588SBarry Smith } 1026d64ed03dSBarry Smith } 10273a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1028d3e5ee88SLois Curfman McInnes Mat tmat; 1029ec8511deSBarry Smith Mat_SeqDense *tmatd; 103087828ca2SBarry Smith PetscScalar *v2; 1031ea709b57SSatish Balay 1032273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1033ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10340de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1035d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 1036*1b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1037d3e5ee88SLois Curfman McInnes } 10386d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10396d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1040d3e5ee88SLois Curfman McInnes *matout = tmat; 104148b35521SBarry Smith } 10423a40ed3dSBarry Smith PetscFunctionReturn(0); 1043289bc588SBarry Smith } 1044289bc588SBarry Smith 10454a2ae208SSatish Balay #undef __FUNCT__ 10464a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 10478f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1048289bc588SBarry Smith { 1049c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1050c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1051*1b807ce4Svictorle int ierr,i,j; 105287828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1053273d9f13SBarry Smith PetscTruth flag; 10549ea5d5aeSSatish Balay 10553a40ed3dSBarry Smith PetscFunctionBegin; 1056273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr); 1057273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 1058273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1059273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1060*1b807ce4Svictorle for (i=0; i<A1->m; i++) { 1061*1b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1062*1b807ce4Svictorle for (j=0; j<A1->n; j++) { 10633a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1064*1b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 1065*1b807ce4Svictorle } 1066289bc588SBarry Smith } 106777c4ece6SBarry Smith *flg = PETSC_TRUE; 10683a40ed3dSBarry Smith PetscFunctionReturn(0); 1069289bc588SBarry Smith } 1070289bc588SBarry Smith 10714a2ae208SSatish Balay #undef __FUNCT__ 10724a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 10738f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1074289bc588SBarry Smith { 1075c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10767a97a34bSBarry Smith int ierr,i,n,len; 107787828ca2SBarry Smith PetscScalar *x,zero = 0.0; 107844cd7ae7SLois Curfman McInnes 10793a40ed3dSBarry Smith PetscFunctionBegin; 10807a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10817a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1082600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1083273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1084273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 108544cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 1086*1b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1087289bc588SBarry Smith } 10887a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10893a40ed3dSBarry Smith PetscFunctionReturn(0); 1090289bc588SBarry Smith } 1091289bc588SBarry Smith 10924a2ae208SSatish Balay #undef __FUNCT__ 10934a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 10948f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1095289bc588SBarry Smith { 1096c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 109787828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1098273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 109955659b69SBarry Smith 11003a40ed3dSBarry Smith PetscFunctionBegin; 110128988994SBarry Smith if (ll) { 11027a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1103600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1104273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1105da3a660dSBarry Smith for (i=0; i<m; i++) { 1106da3a660dSBarry Smith x = l[i]; 1107da3a660dSBarry Smith v = mat->v + i; 1108da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1109da3a660dSBarry Smith } 11107a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1111b0a32e0cSBarry Smith PetscLogFlops(n*m); 1112da3a660dSBarry Smith } 111328988994SBarry Smith if (rr) { 11147a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1115600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1116273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1117da3a660dSBarry Smith for (i=0; i<n; i++) { 1118da3a660dSBarry Smith x = r[i]; 1119da3a660dSBarry Smith v = mat->v + i*m; 1120da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1121da3a660dSBarry Smith } 11227a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1123b0a32e0cSBarry Smith PetscLogFlops(n*m); 1124da3a660dSBarry Smith } 11253a40ed3dSBarry Smith PetscFunctionReturn(0); 1126289bc588SBarry Smith } 1127289bc588SBarry Smith 11284a2ae208SSatish Balay #undef __FUNCT__ 11294a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1130064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1131289bc588SBarry Smith { 1132c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 113387828ca2SBarry Smith PetscScalar *v = mat->v; 1134329f5518SBarry Smith PetscReal sum = 0.0; 1135289bc588SBarry Smith int i,j; 113655659b69SBarry Smith 11373a40ed3dSBarry Smith PetscFunctionBegin; 1138289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1139*1b807ce4Svictorle if (mat->lda>A->m) SETERRQ(1,"Can not handle LDA"); 1140273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1141aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1142329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1143289bc588SBarry Smith #else 1144289bc588SBarry Smith sum += (*v)*(*v); v++; 1145289bc588SBarry Smith #endif 1146289bc588SBarry Smith } 1147064f8208SBarry Smith *nrm = sqrt(sum); 1148b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11493a40ed3dSBarry Smith } else if (type == NORM_1) { 1150064f8208SBarry Smith *nrm = 0.0; 1151273d9f13SBarry Smith for (j=0; j<A->n; j++) { 1152*1b807ce4Svictorle v = mat->v + j*mat->lda; 1153289bc588SBarry Smith sum = 0.0; 1154273d9f13SBarry Smith for (i=0; i<A->m; i++) { 115533a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1156289bc588SBarry Smith } 1157064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1158289bc588SBarry Smith } 1159b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11603a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1161064f8208SBarry Smith *nrm = 0.0; 1162273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1163289bc588SBarry Smith v = mat->v + j; 1164289bc588SBarry Smith sum = 0.0; 1165273d9f13SBarry Smith for (i=0; i<A->n; i++) { 1166*1b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1167289bc588SBarry Smith } 1168064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1169289bc588SBarry Smith } 1170b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11713a40ed3dSBarry Smith } else { 117229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1173289bc588SBarry Smith } 11743a40ed3dSBarry Smith PetscFunctionReturn(0); 1175289bc588SBarry Smith } 1176289bc588SBarry Smith 11774a2ae208SSatish Balay #undef __FUNCT__ 11784a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 11798f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1180289bc588SBarry Smith { 1181c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 118267e560aaSBarry Smith 11833a40ed3dSBarry Smith PetscFunctionBegin; 1184b5a2b587SKris Buschelman switch (op) { 1185b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1186b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1187b5a2b587SKris Buschelman break; 1188b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1189b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1190b5a2b587SKris Buschelman break; 1191b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1192b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1193b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1194b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1195b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1196b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1197b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1198b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1199b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1200b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1201b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1202b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1203b5a2b587SKris Buschelman break; 1204b5a2b587SKris Buschelman default: 120529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 12063a40ed3dSBarry Smith } 12073a40ed3dSBarry Smith PetscFunctionReturn(0); 1208289bc588SBarry Smith } 1209289bc588SBarry Smith 12104a2ae208SSatish Balay #undef __FUNCT__ 12114a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 12128f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 12136f0a148fSBarry Smith { 1214ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1215549d3d68SSatish Balay int ierr; 12163a40ed3dSBarry Smith 12173a40ed3dSBarry Smith PetscFunctionBegin; 1218*1b807ce4Svictorle if (l->lda>A->m) SETERRQ(1,"Can not handle LDA"); 121987828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 12203a40ed3dSBarry Smith PetscFunctionReturn(0); 12216f0a148fSBarry Smith } 12226f0a148fSBarry Smith 12234a2ae208SSatish Balay #undef __FUNCT__ 12244a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 12258f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 12264e220ebcSLois Curfman McInnes { 12273a40ed3dSBarry Smith PetscFunctionBegin; 12284e220ebcSLois Curfman McInnes *bs = 1; 12293a40ed3dSBarry Smith PetscFunctionReturn(0); 12304e220ebcSLois Curfman McInnes } 12314e220ebcSLois Curfman McInnes 12324a2ae208SSatish Balay #undef __FUNCT__ 12334a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 123487828ca2SBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag) 12356f0a148fSBarry Smith { 1236ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1237273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 123887828ca2SBarry Smith PetscScalar *slot; 123955659b69SBarry Smith 12403a40ed3dSBarry Smith PetscFunctionBegin; 1241e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 124278b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12436f0a148fSBarry Smith for (i=0; i<N; i++) { 12446f0a148fSBarry Smith slot = l->v + rows[i]; 12456f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 12466f0a148fSBarry Smith } 12476f0a148fSBarry Smith if (diag) { 12486f0a148fSBarry Smith for (i=0; i<N; i++) { 12496f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1250260925b8SBarry Smith *slot = *diag; 12516f0a148fSBarry Smith } 12526f0a148fSBarry Smith } 1253260925b8SBarry Smith ISRestoreIndices(is,&rows); 12543a40ed3dSBarry Smith PetscFunctionReturn(0); 12556f0a148fSBarry Smith } 1256557bce09SLois Curfman McInnes 12574a2ae208SSatish Balay #undef __FUNCT__ 12584a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 125987828ca2SBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar **array) 126064e87e97SBarry Smith { 1261c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12623a40ed3dSBarry Smith 12633a40ed3dSBarry Smith PetscFunctionBegin; 126464e87e97SBarry Smith *array = mat->v; 12653a40ed3dSBarry Smith PetscFunctionReturn(0); 126664e87e97SBarry Smith } 12670754003eSLois Curfman McInnes 12684a2ae208SSatish Balay #undef __FUNCT__ 12694a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 127087828ca2SBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar **array) 1271ff14e315SSatish Balay { 12723a40ed3dSBarry Smith PetscFunctionBegin; 127309b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12743a40ed3dSBarry Smith PetscFunctionReturn(0); 1275ff14e315SSatish Balay } 12760754003eSLois Curfman McInnes 12774a2ae208SSatish Balay #undef __FUNCT__ 12784a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 12797b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12800754003eSLois Curfman McInnes { 1281c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1282273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 128387828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 12840754003eSLois Curfman McInnes Mat newmat; 12850754003eSLois Curfman McInnes 12863a40ed3dSBarry Smith PetscFunctionBegin; 128778b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 128878b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1289e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1290e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 12910754003eSLois Curfman McInnes 1292182d2002SSatish Balay /* Check submatrixcall */ 1293182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1294182d2002SSatish Balay int n_cols,n_rows; 1295182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 129629bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1297182d2002SSatish Balay newmat = *B; 1298182d2002SSatish Balay } else { 12990754003eSLois Curfman McInnes /* Create and fill new matrix */ 1300b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1301182d2002SSatish Balay } 1302182d2002SSatish Balay 1303182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1304182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1305182d2002SSatish Balay 1306182d2002SSatish Balay for (i=0; i<ncols; i++) { 1307182d2002SSatish Balay av = v + m*icol[i]; 1308182d2002SSatish Balay for (j=0; j<nrows; j++) { 1309182d2002SSatish Balay *bv++ = av[irow[j]]; 13100754003eSLois Curfman McInnes } 13110754003eSLois Curfman McInnes } 1312182d2002SSatish Balay 1313182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13146d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13156d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13160754003eSLois Curfman McInnes 13170754003eSLois Curfman McInnes /* Free work space */ 131878b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 131978b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1320182d2002SSatish Balay *B = newmat; 13213a40ed3dSBarry Smith PetscFunctionReturn(0); 13220754003eSLois Curfman McInnes } 13230754003eSLois Curfman McInnes 13244a2ae208SSatish Balay #undef __FUNCT__ 13254a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 13267b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1327905e6a2fSBarry Smith { 1328905e6a2fSBarry Smith int ierr,i; 1329905e6a2fSBarry Smith 13303a40ed3dSBarry Smith PetscFunctionBegin; 1331905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1332b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1333905e6a2fSBarry Smith } 1334905e6a2fSBarry Smith 1335905e6a2fSBarry Smith for (i=0; i<n; i++) { 13366a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1337905e6a2fSBarry Smith } 13383a40ed3dSBarry Smith PetscFunctionReturn(0); 1339905e6a2fSBarry Smith } 1340905e6a2fSBarry Smith 13414a2ae208SSatish Balay #undef __FUNCT__ 13424a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1343cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 13444b0e389bSBarry Smith { 13454b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 13463a40ed3dSBarry Smith int ierr; 1347273d9f13SBarry Smith PetscTruth flag; 13483a40ed3dSBarry Smith 13493a40ed3dSBarry Smith PetscFunctionBegin; 1350273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr); 1351273d9f13SBarry Smith if (!flag) { 1352cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 13533a40ed3dSBarry Smith PetscFunctionReturn(0); 13543a40ed3dSBarry Smith } 1355273d9f13SBarry Smith if (A->m != B->m || A->n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1356*1b807ce4Svictorle if (a->lda>A->m || b->lda>B->m) SETERRQ(1,"Can not handle LDA"); 135787828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1358273d9f13SBarry Smith PetscFunctionReturn(0); 1359273d9f13SBarry Smith } 1360273d9f13SBarry Smith 13614a2ae208SSatish Balay #undef __FUNCT__ 13624a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1363273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1364273d9f13SBarry Smith { 1365273d9f13SBarry Smith int ierr; 1366273d9f13SBarry Smith 1367273d9f13SBarry Smith PetscFunctionBegin; 1368273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 13693a40ed3dSBarry Smith PetscFunctionReturn(0); 13704b0e389bSBarry Smith } 13714b0e389bSBarry Smith 1372289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1373a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1374905e6a2fSBarry Smith MatGetRow_SeqDense, 1375905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1376905e6a2fSBarry Smith MatMult_SeqDense, 1377905e6a2fSBarry Smith MatMultAdd_SeqDense, 13787c922b88SBarry Smith MatMultTranspose_SeqDense, 13797c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1380905e6a2fSBarry Smith MatSolve_SeqDense, 1381905e6a2fSBarry Smith MatSolveAdd_SeqDense, 13827c922b88SBarry Smith MatSolveTranspose_SeqDense, 13837c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1384905e6a2fSBarry Smith MatLUFactor_SeqDense, 1385905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1386ec8511deSBarry Smith MatRelax_SeqDense, 1387ec8511deSBarry Smith MatTranspose_SeqDense, 1388905e6a2fSBarry Smith MatGetInfo_SeqDense, 1389905e6a2fSBarry Smith MatEqual_SeqDense, 1390905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1391905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1392905e6a2fSBarry Smith MatNorm_SeqDense, 1393905e6a2fSBarry Smith 0, 1394905e6a2fSBarry Smith 0, 1395905e6a2fSBarry Smith 0, 1396905e6a2fSBarry Smith MatSetOption_SeqDense, 1397905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1398905e6a2fSBarry Smith MatZeroRows_SeqDense, 1399905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1400905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1401905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1402905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1403273d9f13SBarry Smith MatSetUpPreallocation_SeqDense, 1404273d9f13SBarry Smith 0, 1405905e6a2fSBarry Smith 0, 1406905e6a2fSBarry Smith MatGetArray_SeqDense, 1407905e6a2fSBarry Smith MatRestoreArray_SeqDense, 14085609ef8eSBarry Smith MatDuplicate_SeqDense, 1409a5ae1ecdSBarry Smith 0, 1410a5ae1ecdSBarry Smith 0, 1411a5ae1ecdSBarry Smith 0, 1412a5ae1ecdSBarry Smith 0, 1413a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1414a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1415a5ae1ecdSBarry Smith 0, 14164b0e389bSBarry Smith MatGetValues_SeqDense, 1417a5ae1ecdSBarry Smith MatCopy_SeqDense, 1418a5ae1ecdSBarry Smith 0, 1419a5ae1ecdSBarry Smith MatScale_SeqDense, 1420a5ae1ecdSBarry Smith 0, 1421a5ae1ecdSBarry Smith 0, 1422a5ae1ecdSBarry Smith 0, 1423a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1424a5ae1ecdSBarry Smith 0, 1425a5ae1ecdSBarry Smith 0, 1426a5ae1ecdSBarry Smith 0, 1427a5ae1ecdSBarry Smith 0, 1428a5ae1ecdSBarry Smith 0, 1429a5ae1ecdSBarry Smith 0, 1430a5ae1ecdSBarry Smith 0, 1431a5ae1ecdSBarry Smith 0, 1432a5ae1ecdSBarry Smith 0, 1433a5ae1ecdSBarry Smith 0, 1434e03a110bSBarry Smith MatDestroy_SeqDense, 1435e03a110bSBarry Smith MatView_SeqDense, 14368a124369SBarry Smith MatGetPetscMaps_Petsc}; 143790ace30eSBarry Smith 14384a2ae208SSatish Balay #undef __FUNCT__ 14394a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 14404b828684SBarry Smith /*@C 1441fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1442d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1443d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1444289bc588SBarry Smith 1445db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1446db81eaa0SLois Curfman McInnes 144720563c6bSBarry Smith Input Parameters: 1448db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 14490c775827SLois Curfman McInnes . m - number of rows 145018f449edSLois Curfman McInnes . n - number of columns 1451db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1452dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 145320563c6bSBarry Smith 145420563c6bSBarry Smith Output Parameter: 145544cd7ae7SLois Curfman McInnes . A - the matrix 145620563c6bSBarry Smith 1457b259b22eSLois Curfman McInnes Notes: 145818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 145918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1460b4fd4287SBarry Smith set data=PETSC_NULL. 146118f449edSLois Curfman McInnes 1462027ccd11SLois Curfman McInnes Level: intermediate 1463027ccd11SLois Curfman McInnes 1464dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1465d65003e9SLois Curfman McInnes 1466db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 146720563c6bSBarry Smith @*/ 146887828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1469289bc588SBarry Smith { 1470273d9f13SBarry Smith int ierr; 14713b2fbd54SBarry Smith 14723a40ed3dSBarry Smith PetscFunctionBegin; 1473273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1474273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1475273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1476273d9f13SBarry Smith PetscFunctionReturn(0); 1477273d9f13SBarry Smith } 1478273d9f13SBarry Smith 14794a2ae208SSatish Balay #undef __FUNCT__ 14804a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1481273d9f13SBarry Smith /*@C 1482273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1483273d9f13SBarry Smith 1484273d9f13SBarry Smith Collective on MPI_Comm 1485273d9f13SBarry Smith 1486273d9f13SBarry Smith Input Parameters: 1487273d9f13SBarry Smith + A - the matrix 1488273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1489273d9f13SBarry Smith 1490273d9f13SBarry Smith Notes: 1491273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1492273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1493273d9f13SBarry Smith set data=PETSC_NULL. 1494273d9f13SBarry Smith 1495273d9f13SBarry Smith Level: intermediate 1496273d9f13SBarry Smith 1497273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1498273d9f13SBarry Smith 1499273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1500273d9f13SBarry Smith @*/ 150187828ca2SBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data) 1502273d9f13SBarry Smith { 1503273d9f13SBarry Smith Mat_SeqDense *b; 1504273d9f13SBarry Smith int ierr; 1505273d9f13SBarry Smith PetscTruth flg2; 1506273d9f13SBarry Smith 1507273d9f13SBarry Smith PetscFunctionBegin; 1508273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flg2);CHKERRQ(ierr); 1509273d9f13SBarry Smith if (!flg2) PetscFunctionReturn(0); 1510273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1511273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1512273d9f13SBarry Smith if (!data) { 151387828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 151487828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1515273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 151687828ca2SBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1517273d9f13SBarry Smith } else { /* user-allocated storage */ 1518273d9f13SBarry Smith b->v = data; 1519273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1520273d9f13SBarry Smith } 1521273d9f13SBarry Smith PetscFunctionReturn(0); 1522273d9f13SBarry Smith } 1523273d9f13SBarry Smith 1524*1b807ce4Svictorle #undef __FUNCT__ 1525*1b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 1526*1b807ce4Svictorle /*@C 1527*1b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 1528*1b807ce4Svictorle 1529*1b807ce4Svictorle Input parameter: 1530*1b807ce4Svictorle + A - the matrix 1531*1b807ce4Svictorle - lda - the leading dimension 1532*1b807ce4Svictorle 1533*1b807ce4Svictorle Notes: 1534*1b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 1535*1b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 1536*1b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 1537*1b807ce4Svictorle 1538*1b807ce4Svictorle Level: intermediate 1539*1b807ce4Svictorle 1540*1b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 1541*1b807ce4Svictorle 1542*1b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 1543*1b807ce4Svictorle @*/ 1544*1b807ce4Svictorle int MatSeqDenseSetLDA(Mat B,int lda) 1545*1b807ce4Svictorle { 1546*1b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1547*1b807ce4Svictorle PetscFunctionBegin; 1548*1b807ce4Svictorle if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension"); 1549*1b807ce4Svictorle b->lda = lda; 1550*1b807ce4Svictorle PetscFunctionReturn(0); 1551*1b807ce4Svictorle } 1552*1b807ce4Svictorle 1553273d9f13SBarry Smith EXTERN_C_BEGIN 15544a2ae208SSatish Balay #undef __FUNCT__ 15554a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1556273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1557273d9f13SBarry Smith { 1558273d9f13SBarry Smith Mat_SeqDense *b; 1559273d9f13SBarry Smith int ierr,size; 1560273d9f13SBarry Smith 1561273d9f13SBarry Smith PetscFunctionBegin; 1562273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 156329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 156455659b69SBarry Smith 1565273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1566273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1567273d9f13SBarry Smith 1568b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1569549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1570549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 157144cd7ae7SLois Curfman McInnes B->factor = 0; 157290f02eecSBarry Smith B->mapping = 0; 1573b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 157444cd7ae7SLois Curfman McInnes B->data = (void*)b; 157518f449edSLois Curfman McInnes 15768a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 15778a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1578a5ae1ecdSBarry Smith 157944cd7ae7SLois Curfman McInnes b->pivots = 0; 1580273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1581273d9f13SBarry Smith b->v = 0; 1582*1b807ce4Svictorle b->lda = B->m; 15834e220ebcSLois Curfman McInnes 15843a40ed3dSBarry Smith PetscFunctionReturn(0); 1585289bc588SBarry Smith } 1586273d9f13SBarry Smith EXTERN_C_END 1587