167e560aaSBarry Smith /* 267e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 367e560aaSBarry Smith */ 4289bc588SBarry Smith 570f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 6b4c862e0SSatish Balay #include "petscblaslapack.h" 7289bc588SBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 10dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str) 111987afe7SBarry Smith { 121987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 13a5ce6ee0Svictorle int N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1; 143a40ed3dSBarry Smith 153a40ed3dSBarry Smith PetscFunctionBegin; 16a5ce6ee0Svictorle if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 17a5ce6ee0Svictorle if (ldax>m || lday>m) { 18a5ce6ee0Svictorle for (j=0; j<X->n; j++) { 19268466fbSBarry Smith BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one); 20a5ce6ee0Svictorle } 21a5ce6ee0Svictorle } else { 22268466fbSBarry Smith BLaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one); 23a5ce6ee0Svictorle } 24b0a32e0cSBarry Smith PetscLogFlops(2*N-1); 253a40ed3dSBarry Smith PetscFunctionReturn(0); 261987afe7SBarry Smith } 271987afe7SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 30dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 31289bc588SBarry Smith { 324e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 33273d9f13SBarry Smith int i,N = A->m*A->n,count = 0; 3487828ca2SBarry Smith PetscScalar *v = a->v; 353a40ed3dSBarry Smith 363a40ed3dSBarry Smith PetscFunctionBegin; 37289bc588SBarry Smith for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 384e220ebcSLois Curfman McInnes 39273d9f13SBarry Smith info->rows_global = (double)A->m; 40273d9f13SBarry Smith info->columns_global = (double)A->n; 41273d9f13SBarry Smith info->rows_local = (double)A->m; 42273d9f13SBarry Smith info->columns_local = (double)A->n; 434e220ebcSLois Curfman McInnes info->block_size = 1.0; 444e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 454e220ebcSLois Curfman McInnes info->nz_used = (double)count; 464e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 474e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 484e220ebcSLois Curfman McInnes info->mallocs = 0; 494e220ebcSLois Curfman McInnes info->memory = A->mem; 504e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 514e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 524e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 534e220ebcSLois Curfman McInnes 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55289bc588SBarry Smith } 56289bc588SBarry Smith 574a2ae208SSatish Balay #undef __FUNCT__ 584a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 59dfbe8321SBarry Smith PetscErrorCode MatScale_SeqDense(const PetscScalar *alpha,Mat A) 6080cd9d93SLois Curfman McInnes { 61273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 62a5ce6ee0Svictorle int one = 1,lda = a->lda,j,nz; 6380cd9d93SLois Curfman McInnes 643a40ed3dSBarry Smith PetscFunctionBegin; 65a5ce6ee0Svictorle if (lda>A->m) { 66a5ce6ee0Svictorle nz = A->m; 67a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 68268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one); 69a5ce6ee0Svictorle } 70a5ce6ee0Svictorle } else { 71273d9f13SBarry Smith nz = A->m*A->n; 72268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 73a5ce6ee0Svictorle } 74b0a32e0cSBarry Smith PetscLogFlops(nz); 753a40ed3dSBarry Smith PetscFunctionReturn(0); 7680cd9d93SLois Curfman McInnes } 7780cd9d93SLois Curfman McInnes 78289bc588SBarry Smith /* ---------------------------------------------------------------*/ 796831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 806831982aSBarry Smith rather than put it in the Mat->row slot.*/ 814a2ae208SSatish Balay #undef __FUNCT__ 824a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 83dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 84289bc588SBarry Smith { 85a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 8681f479a6Svictorle PetscFunctionBegin; 87a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 88a07ab388SBarry Smith #else 89c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 90*6849ba73SBarry Smith PetscErrorCode ierr; 91*6849ba73SBarry Smith int info; 9255659b69SBarry Smith 933a40ed3dSBarry Smith PetscFunctionBegin; 94289bc588SBarry Smith if (!mat->pivots) { 95b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr); 96b0a32e0cSBarry Smith PetscLogObjectMemory(A,A->m*sizeof(int)); 97289bc588SBarry Smith } 98f1af5d2fSBarry Smith A->factor = FACTOR_LU; 99273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 1001b807ce4Svictorle LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info); 10129bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 10229bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 103b0a32e0cSBarry Smith PetscLogFlops((2*A->n*A->n*A->n)/3); 104a07ab388SBarry Smith #endif 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 106289bc588SBarry Smith } 1076ee01492SSatish Balay 1084a2ae208SSatish Balay #undef __FUNCT__ 1094a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 110dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11102cad45dSBarry Smith { 11202cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 113*6849ba73SBarry Smith PetscErrorCode ierr; 114*6849ba73SBarry Smith int lda = mat->lda,j,m; 11502cad45dSBarry Smith Mat newi; 11602cad45dSBarry Smith 1173a40ed3dSBarry Smith PetscFunctionBegin; 1185c5985e7SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&newi);CHKERRQ(ierr); 1195c5985e7SKris Buschelman ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr); 1205c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 1215609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 122a5ce6ee0Svictorle l = (Mat_SeqDense*)newi->data; 123a5ce6ee0Svictorle if (lda>A->m) { 124a5ce6ee0Svictorle m = A->m; 125a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 126a5ce6ee0Svictorle ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 127a5ce6ee0Svictorle } 128a5ce6ee0Svictorle } else { 12987828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 13002cad45dSBarry Smith } 131a5ce6ee0Svictorle } 13202cad45dSBarry Smith newi->assembled = PETSC_TRUE; 13302cad45dSBarry Smith *newmat = newi; 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 13502cad45dSBarry Smith } 13602cad45dSBarry Smith 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 139dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 140289bc588SBarry Smith { 141dfbe8321SBarry Smith PetscErrorCode ierr; 1423a40ed3dSBarry Smith 1433a40ed3dSBarry Smith PetscFunctionBegin; 1445609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1453a40ed3dSBarry Smith PetscFunctionReturn(0); 146289bc588SBarry Smith } 1476ee01492SSatish Balay 1484a2ae208SSatish Balay #undef __FUNCT__ 1494a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 150dfbe8321SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 151289bc588SBarry Smith { 15202cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 153*6849ba73SBarry Smith PetscErrorCode ierr; 154*6849ba73SBarry Smith int lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j; 1554482741eSBarry Smith MatFactorInfo info; 1563a40ed3dSBarry Smith 1573a40ed3dSBarry Smith PetscFunctionBegin; 15802cad45dSBarry Smith /* copy the numerical values */ 1591b807ce4Svictorle if (lda1>m || lda2>m ) { 1601b807ce4Svictorle for (j=0; j<n; j++) { 1611b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1621b807ce4Svictorle } 1631b807ce4Svictorle } else { 16487828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1651b807ce4Svictorle } 16602cad45dSBarry Smith (*fact)->factor = 0; 1674482741eSBarry Smith ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr); 1683a40ed3dSBarry Smith PetscFunctionReturn(0); 169289bc588SBarry Smith } 1706ee01492SSatish Balay 1714a2ae208SSatish Balay #undef __FUNCT__ 1724a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 173dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 174289bc588SBarry Smith { 175dfbe8321SBarry Smith PetscErrorCode ierr; 1763a40ed3dSBarry Smith 1773a40ed3dSBarry Smith PetscFunctionBegin; 1783a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1793a40ed3dSBarry Smith PetscFunctionReturn(0); 180289bc588SBarry Smith } 1816ee01492SSatish Balay 1824a2ae208SSatish Balay #undef __FUNCT__ 1834a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 184dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 185289bc588SBarry Smith { 186a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 187a07ab388SBarry Smith PetscFunctionBegin; 188a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 189a07ab388SBarry Smith #else 190c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 191*6849ba73SBarry Smith PetscErrorCode ierr; 192*6849ba73SBarry Smith int info; 19355659b69SBarry Smith 1943a40ed3dSBarry Smith PetscFunctionBegin; 195752f5567SLois Curfman McInnes if (mat->pivots) { 196606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 197b0a32e0cSBarry Smith PetscLogObjectMemory(A,-A->m*sizeof(int)); 198752f5567SLois Curfman McInnes mat->pivots = 0; 199752f5567SLois Curfman McInnes } 200273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2011b807ce4Svictorle LApotrf_("L",&A->n,mat->v,&mat->lda,&info); 20229bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 203c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 204b0a32e0cSBarry Smith PetscLogFlops((A->n*A->n*A->n)/3); 205a07ab388SBarry Smith #endif 2063a40ed3dSBarry Smith PetscFunctionReturn(0); 207289bc588SBarry Smith } 208289bc588SBarry Smith 2094a2ae208SSatish Balay #undef __FUNCT__ 2100b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 211dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 2120b4b3355SBarry Smith { 213dfbe8321SBarry Smith PetscErrorCode ierr; 21415e8a5b3SHong Zhang MatFactorInfo info; 2150b4b3355SBarry Smith 2160b4b3355SBarry Smith PetscFunctionBegin; 21715e8a5b3SHong Zhang info.fill = 1.0; 21815e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 2190b4b3355SBarry Smith PetscFunctionReturn(0); 2200b4b3355SBarry Smith } 2210b4b3355SBarry Smith 2220b4b3355SBarry Smith #undef __FUNCT__ 2234a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 224dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 225289bc588SBarry Smith { 226c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 227*6849ba73SBarry Smith PetscErrorCode ierr; 228*6849ba73SBarry Smith int one = 1,info; 22987828ca2SBarry Smith PetscScalar *x,*y; 23067e560aaSBarry Smith 2313a40ed3dSBarry Smith PetscFunctionBegin; 232273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2331ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2341ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 23587828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 236c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 237ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 23880444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 239ae7cfcebSSatish Balay #else 2401b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2412ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 242ae7cfcebSSatish Balay #endif 2437a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 244ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 24580444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 246ae7cfcebSSatish Balay #else 2471b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2482ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 249ae7cfcebSSatish Balay #endif 250289bc588SBarry Smith } 25129bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2521ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2531ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 254b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2553a40ed3dSBarry Smith PetscFunctionReturn(0); 256289bc588SBarry Smith } 2576ee01492SSatish Balay 2584a2ae208SSatish Balay #undef __FUNCT__ 2594a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 260dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 261da3a660dSBarry Smith { 262c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 263dfbe8321SBarry Smith PetscErrorCode ierr; 264dfbe8321SBarry Smith int one = 1,info; 26587828ca2SBarry Smith PetscScalar *x,*y; 26667e560aaSBarry Smith 2673a40ed3dSBarry Smith PetscFunctionBegin; 268273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2691ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2701ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 27187828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 272752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 273da3a660dSBarry Smith if (mat->pivots) { 274ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 27580444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 276ae7cfcebSSatish Balay #else 2771b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2782ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 279ae7cfcebSSatish Balay #endif 2807a97a34bSBarry Smith } else { 281ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 28280444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 283ae7cfcebSSatish Balay #else 2841b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2852ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 286ae7cfcebSSatish Balay #endif 287da3a660dSBarry Smith } 2881ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2891ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 290b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2913a40ed3dSBarry Smith PetscFunctionReturn(0); 292da3a660dSBarry Smith } 2936ee01492SSatish Balay 2944a2ae208SSatish Balay #undef __FUNCT__ 2954a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 296dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 297da3a660dSBarry Smith { 298c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 299dfbe8321SBarry Smith PetscErrorCode ierr; 300dfbe8321SBarry Smith int one = 1,info; 30187828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 302da3a660dSBarry Smith Vec tmp = 0; 30367e560aaSBarry Smith 3043a40ed3dSBarry Smith PetscFunctionBegin; 3051ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3061ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 307273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 308da3a660dSBarry Smith if (yy == zz) { 30978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 310b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 31178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 312da3a660dSBarry Smith } 31387828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 314752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 315da3a660dSBarry Smith if (mat->pivots) { 316ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 31780444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 318ae7cfcebSSatish Balay #else 3191b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 3202ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 321ae7cfcebSSatish Balay #endif 322a8c6a408SBarry Smith } else { 323ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 32480444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 325ae7cfcebSSatish Balay #else 3261b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 3272ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 328ae7cfcebSSatish Balay #endif 329da3a660dSBarry Smith } 330a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 331a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 3321ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3331ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 334b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 336da3a660dSBarry Smith } 33767e560aaSBarry Smith 3384a2ae208SSatish Balay #undef __FUNCT__ 3394a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 340dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 341da3a660dSBarry Smith { 342c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 343*6849ba73SBarry Smith PetscErrorCode ierr; 344*6849ba73SBarry Smith int one = 1,info; 34587828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 346da3a660dSBarry Smith Vec tmp; 34767e560aaSBarry Smith 3483a40ed3dSBarry Smith PetscFunctionBegin; 349273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3501ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3511ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 352da3a660dSBarry Smith if (yy == zz) { 35378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 354b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 35578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 356da3a660dSBarry Smith } 35787828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 358752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 359da3a660dSBarry Smith if (mat->pivots) { 360ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 36180444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 362ae7cfcebSSatish Balay #else 3631b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 3642ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 365ae7cfcebSSatish Balay #endif 3663a40ed3dSBarry Smith } else { 367ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 36880444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 369ae7cfcebSSatish Balay #else 3701b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 3712ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 372ae7cfcebSSatish Balay #endif 373da3a660dSBarry Smith } 37490f02eecSBarry Smith if (tmp) { 37590f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 37690f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3773a40ed3dSBarry Smith } else { 37890f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 37990f02eecSBarry Smith } 3801ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3811ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 382b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3833a40ed3dSBarry Smith PetscFunctionReturn(0); 384da3a660dSBarry Smith } 385289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3864a2ae208SSatish Balay #undef __FUNCT__ 3874a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 388dfbe8321SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 389c14dc6b6SHong Zhang PetscReal shift,int its,int lits,Vec xx) 390289bc588SBarry Smith { 391c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 39287828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 393dfbe8321SBarry Smith PetscErrorCode ierr; 394dfbe8321SBarry Smith int m = A->m,i; 395aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 396bc1b551cSSatish Balay int o = 1; 397bc1b551cSSatish Balay #endif 398289bc588SBarry Smith 3993a40ed3dSBarry Smith PetscFunctionBegin; 400289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 4013a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 402bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 403289bc588SBarry Smith } 4041ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4051ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 406b965ef7fSBarry Smith its = its*lits; 40791723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 408289bc588SBarry Smith while (its--) { 409289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 410289bc588SBarry Smith for (i=0; i<m; i++) { 411aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 412f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 413f1747703SBarry Smith not happy about returning a double complex */ 414f1747703SBarry Smith int _i; 41587828ca2SBarry Smith PetscScalar sum = b[i]; 416f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4173f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 418f1747703SBarry Smith } 419f1747703SBarry Smith xt = sum; 420f1747703SBarry Smith #else 421289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 422f1747703SBarry Smith #endif 42355a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 424289bc588SBarry Smith } 425289bc588SBarry Smith } 426289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 427289bc588SBarry Smith for (i=m-1; i>=0; i--) { 428aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 429f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 430f1747703SBarry Smith not happy about returning a double complex */ 431f1747703SBarry Smith int _i; 43287828ca2SBarry Smith PetscScalar sum = b[i]; 433f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4343f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 435f1747703SBarry Smith } 436f1747703SBarry Smith xt = sum; 437f1747703SBarry Smith #else 438289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 439f1747703SBarry Smith #endif 44055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 441289bc588SBarry Smith } 442289bc588SBarry Smith } 443289bc588SBarry Smith } 4441ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4451ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4463a40ed3dSBarry Smith PetscFunctionReturn(0); 447289bc588SBarry Smith } 448289bc588SBarry Smith 449289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4504a2ae208SSatish Balay #undef __FUNCT__ 4514a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 452dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 453289bc588SBarry Smith { 454c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 45587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 456dfbe8321SBarry Smith PetscErrorCode ierr; 457dfbe8321SBarry Smith int _One=1; 458ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4593a40ed3dSBarry Smith 4603a40ed3dSBarry Smith PetscFunctionBegin; 461273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4621ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4631ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4641b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4651ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4661ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 467b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->n); 4683a40ed3dSBarry Smith PetscFunctionReturn(0); 469289bc588SBarry Smith } 4706ee01492SSatish Balay 4714a2ae208SSatish Balay #undef __FUNCT__ 4724a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 473dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 474289bc588SBarry Smith { 475c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 47687828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 477dfbe8321SBarry Smith PetscErrorCode ierr; 478dfbe8321SBarry Smith int _One=1; 4793a40ed3dSBarry Smith 4803a40ed3dSBarry Smith PetscFunctionBegin; 481273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4821ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4831ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4841b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4851ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4861ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 487b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->m); 4883a40ed3dSBarry Smith PetscFunctionReturn(0); 489289bc588SBarry Smith } 4906ee01492SSatish Balay 4914a2ae208SSatish Balay #undef __FUNCT__ 4924a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 493dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 494289bc588SBarry Smith { 495c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 49687828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 497dfbe8321SBarry Smith PetscErrorCode ierr; 498dfbe8321SBarry Smith int _One=1; 4993a40ed3dSBarry Smith 5003a40ed3dSBarry Smith PetscFunctionBegin; 501273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 502600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5031ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5041ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5051b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5061ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5071ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 508b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 5093a40ed3dSBarry Smith PetscFunctionReturn(0); 510289bc588SBarry Smith } 5116ee01492SSatish Balay 5124a2ae208SSatish Balay #undef __FUNCT__ 5134a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 514dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 515289bc588SBarry Smith { 516c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 51787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 518dfbe8321SBarry Smith PetscErrorCode ierr; 519dfbe8321SBarry Smith int _One=1; 52087828ca2SBarry Smith PetscScalar _DOne=1.0; 5213a40ed3dSBarry Smith 5223a40ed3dSBarry Smith PetscFunctionBegin; 523273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 524600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5251ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5261ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5271b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5281ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5291ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 530b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 5313a40ed3dSBarry Smith PetscFunctionReturn(0); 532289bc588SBarry Smith } 533289bc588SBarry Smith 534289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5354a2ae208SSatish Balay #undef __FUNCT__ 5364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 537dfbe8321SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 538289bc588SBarry Smith { 539c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54087828ca2SBarry Smith PetscScalar *v; 541*6849ba73SBarry Smith PetscErrorCode ierr; 542*6849ba73SBarry Smith int i; 54367e560aaSBarry Smith 5443a40ed3dSBarry Smith PetscFunctionBegin; 545273d9f13SBarry Smith *ncols = A->n; 546289bc588SBarry Smith if (cols) { 547b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 548273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 549289bc588SBarry Smith } 550289bc588SBarry Smith if (vals) { 55187828ca2SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 552289bc588SBarry Smith v = mat->v + row; 5531b807ce4Svictorle for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;} 554289bc588SBarry Smith } 5553a40ed3dSBarry Smith PetscFunctionReturn(0); 556289bc588SBarry Smith } 5576ee01492SSatish Balay 5584a2ae208SSatish Balay #undef __FUNCT__ 5594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 560dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 561289bc588SBarry Smith { 562dfbe8321SBarry Smith PetscErrorCode ierr; 563606d414cSSatish Balay PetscFunctionBegin; 564606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 565606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5663a40ed3dSBarry Smith PetscFunctionReturn(0); 567289bc588SBarry Smith } 568289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5694a2ae208SSatish Balay #undef __FUNCT__ 5704a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 571dfbe8321SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv) 572289bc588SBarry Smith { 573c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 574cddbea37SSatish Balay int i,j,idx=0; 575d6dfbf8fSBarry Smith 5763a40ed3dSBarry Smith PetscFunctionBegin; 577289bc588SBarry Smith if (!mat->roworiented) { 578dbb450caSBarry Smith if (addv == INSERT_VALUES) { 579289bc588SBarry Smith for (j=0; j<n; j++) { 580cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 581aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 582590ac198SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 58358804f6dSBarry Smith #endif 584289bc588SBarry Smith for (i=0; i<m; i++) { 585cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 586aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 5875c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 58858804f6dSBarry Smith #endif 589cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 590289bc588SBarry Smith } 591289bc588SBarry Smith } 5923a40ed3dSBarry Smith } else { 593289bc588SBarry Smith for (j=0; j<n; j++) { 594cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 595aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 596590ac198SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 59758804f6dSBarry Smith #endif 598289bc588SBarry Smith for (i=0; i<m; i++) { 599cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 600aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6015c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 60258804f6dSBarry Smith #endif 603cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 604289bc588SBarry Smith } 605289bc588SBarry Smith } 606289bc588SBarry Smith } 6073a40ed3dSBarry Smith } else { 608dbb450caSBarry Smith if (addv == INSERT_VALUES) { 609e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 610cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 611aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6125c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 61358804f6dSBarry Smith #endif 614e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 615cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 616aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6175c8e7597SSatish Balay if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 61858804f6dSBarry Smith #endif 619cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 620e8d4e0b9SBarry Smith } 621e8d4e0b9SBarry Smith } 6223a40ed3dSBarry Smith } else { 623289bc588SBarry Smith for (i=0; i<m; i++) { 624cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 625aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6265c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 62758804f6dSBarry Smith #endif 628289bc588SBarry Smith for (j=0; j<n; j++) { 629cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 630aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6315c8e7597SSatish Balay if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 63258804f6dSBarry Smith #endif 633cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 634289bc588SBarry Smith } 635289bc588SBarry Smith } 636289bc588SBarry Smith } 637e8d4e0b9SBarry Smith } 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 639289bc588SBarry Smith } 640e8d4e0b9SBarry Smith 6414a2ae208SSatish Balay #undef __FUNCT__ 6424a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 643dfbe8321SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[]) 644ae80bb75SLois Curfman McInnes { 645ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 646ae80bb75SLois Curfman McInnes int i,j; 64787828ca2SBarry Smith PetscScalar *vpt = v; 648ae80bb75SLois Curfman McInnes 6493a40ed3dSBarry Smith PetscFunctionBegin; 650ae80bb75SLois Curfman McInnes /* row-oriented output */ 651ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 652ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 6531b807ce4Svictorle *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 654ae80bb75SLois Curfman McInnes } 655ae80bb75SLois Curfman McInnes } 6563a40ed3dSBarry Smith PetscFunctionReturn(0); 657ae80bb75SLois Curfman McInnes } 658ae80bb75SLois Curfman McInnes 659289bc588SBarry Smith /* -----------------------------------------------------------------*/ 660289bc588SBarry Smith 661e090d566SSatish Balay #include "petscsys.h" 66288e32edaSLois Curfman McInnes 6634a2ae208SSatish Balay #undef __FUNCT__ 6644a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 665dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A) 66688e32edaSLois Curfman McInnes { 66788e32edaSLois Curfman McInnes Mat_SeqDense *a; 66888e32edaSLois Curfman McInnes Mat B; 669*6849ba73SBarry Smith PetscErrorCode ierr; 670*6849ba73SBarry Smith int *scols,i,j,nz,fd,header[4],size; 67188e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 67287828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 67319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 67488e32edaSLois Curfman McInnes 6753a40ed3dSBarry Smith PetscFunctionBegin; 676d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 67729bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 678b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6790752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 680552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 68188e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 68288e32edaSLois Curfman McInnes 68390ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 6845c5985e7SKris Buschelman ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr); 685be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 6865c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 68790ace30eSBarry Smith B = *A; 68890ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6894a41a4d3SSatish Balay v = a->v; 6904a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6914a41a4d3SSatish Balay from row major to column major */ 692b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 69390ace30eSBarry Smith /* read in nonzero values */ 6944a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6954a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6964a41a4d3SSatish Balay for (j=0; j<N; j++) { 6974a41a4d3SSatish Balay for (i=0; i<M; i++) { 6984a41a4d3SSatish Balay *v++ =w[i*N+j]; 6994a41a4d3SSatish Balay } 7004a41a4d3SSatish Balay } 701606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 7026d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7036d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70490ace30eSBarry Smith } else { 70588e32edaSLois Curfman McInnes /* read row lengths */ 706b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 7070752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 70888e32edaSLois Curfman McInnes 70988e32edaSLois Curfman McInnes /* create our matrix */ 7105c5985e7SKris Buschelman ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr); 711be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7125c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 71388e32edaSLois Curfman McInnes B = *A; 71488e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 71588e32edaSLois Curfman McInnes v = a->v; 71688e32edaSLois Curfman McInnes 71788e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 718b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 719b0a32e0cSBarry Smith cols = scols; 7200752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 72187828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 722b0a32e0cSBarry Smith vals = svals; 7230752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 72488e32edaSLois Curfman McInnes 72588e32edaSLois Curfman McInnes /* insert into matrix */ 72688e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 72788e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 72888e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 72988e32edaSLois Curfman McInnes } 730606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 731606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 732606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 73388e32edaSLois Curfman McInnes 7346d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7356d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 73690ace30eSBarry Smith } 7373a40ed3dSBarry Smith PetscFunctionReturn(0); 73888e32edaSLois Curfman McInnes } 73988e32edaSLois Curfman McInnes 740e090d566SSatish Balay #include "petscsys.h" 741932b0c3eSLois Curfman McInnes 7424a2ae208SSatish Balay #undef __FUNCT__ 7434a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 744*6849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 745289bc588SBarry Smith { 746932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 747dfbe8321SBarry Smith PetscErrorCode ierr; 748dfbe8321SBarry Smith int i,j; 749fb9695e5SSatish Balay char *name; 75087828ca2SBarry Smith PetscScalar *v; 751f3ef73ceSBarry Smith PetscViewerFormat format; 752932b0c3eSLois Curfman McInnes 7533a40ed3dSBarry Smith PetscFunctionBegin; 754435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 755b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 756456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7573a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 758fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 759b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 760273d9f13SBarry Smith for (i=0; i<A->m; i++) { 76144cd7ae7SLois Curfman McInnes v = a->v + i; 762b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 763273d9f13SBarry Smith for (j=0; j<A->n; j++) { 764aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 765329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 76662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 767329f5518SBarry Smith } else if (PetscRealPart(*v)) { 76862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7696831982aSBarry Smith } 77080cd9d93SLois Curfman McInnes #else 7716831982aSBarry Smith if (*v) { 77262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 7736831982aSBarry Smith } 77480cd9d93SLois Curfman McInnes #endif 7751b807ce4Svictorle v += a->lda; 77680cd9d93SLois Curfman McInnes } 777b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 77880cd9d93SLois Curfman McInnes } 779b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7803a40ed3dSBarry Smith } else { 781b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 782aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 783ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 78447989497SBarry Smith /* determine if matrix has all real values */ 78547989497SBarry Smith v = a->v; 786201f5f94SSatish Balay for (i=0; i<A->m*A->n; i++) { 787ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 78847989497SBarry Smith } 78947989497SBarry Smith #endif 790fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7913a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 792b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 793fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 794fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 795ffac6cdbSBarry Smith } 796ffac6cdbSBarry Smith 797273d9f13SBarry Smith for (i=0; i<A->m; i++) { 798932b0c3eSLois Curfman McInnes v = a->v + i; 799273d9f13SBarry Smith for (j=0; j<A->n; j++) { 800aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 80147989497SBarry Smith if (allreal) { 8021b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 80347989497SBarry Smith } else { 8041b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 80547989497SBarry Smith } 806289bc588SBarry Smith #else 8071b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 808289bc588SBarry Smith #endif 8091b807ce4Svictorle v += a->lda; 810289bc588SBarry Smith } 811b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 812289bc588SBarry Smith } 813fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 814b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 815ffac6cdbSBarry Smith } 816b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 817da3a660dSBarry Smith } 818b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8193a40ed3dSBarry Smith PetscFunctionReturn(0); 820289bc588SBarry Smith } 821289bc588SBarry Smith 8224a2ae208SSatish Balay #undef __FUNCT__ 8234a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 824*6849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 825932b0c3eSLois Curfman McInnes { 826932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 827*6849ba73SBarry Smith PetscErrorCode ierr; 828*6849ba73SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,nz = m*n; 82987828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 830f3ef73ceSBarry Smith PetscViewerFormat format; 831932b0c3eSLois Curfman McInnes 8323a40ed3dSBarry Smith PetscFunctionBegin; 833b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 83490ace30eSBarry Smith 835b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 836fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 83790ace30eSBarry Smith /* store the matrix as a dense matrix */ 838b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 839552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 84090ace30eSBarry Smith col_lens[1] = m; 84190ace30eSBarry Smith col_lens[2] = n; 84290ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8430752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 844606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 84590ace30eSBarry Smith 84690ace30eSBarry Smith /* write out matrix, by rows */ 84787828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 84890ace30eSBarry Smith v = a->v; 84990ace30eSBarry Smith for (i=0; i<m; i++) { 85090ace30eSBarry Smith for (j=0; j<n; j++) { 85190ace30eSBarry Smith vals[i + j*m] = *v++; 85290ace30eSBarry Smith } 85390ace30eSBarry Smith } 8540752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 855606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 85690ace30eSBarry Smith } else { 857b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 858552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 859932b0c3eSLois Curfman McInnes col_lens[1] = m; 860932b0c3eSLois Curfman McInnes col_lens[2] = n; 861932b0c3eSLois Curfman McInnes col_lens[3] = nz; 862932b0c3eSLois Curfman McInnes 863932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 864932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8650752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 866932b0c3eSLois Curfman McInnes 867932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 868932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 869932b0c3eSLois Curfman McInnes ict = 0; 870932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 871932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 872932b0c3eSLois Curfman McInnes } 8730752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 874606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 875932b0c3eSLois Curfman McInnes 876932b0c3eSLois Curfman McInnes /* store nonzero values */ 87787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 878932b0c3eSLois Curfman McInnes ict = 0; 879932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 880932b0c3eSLois Curfman McInnes v = a->v + i; 881932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 8821b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 883932b0c3eSLois Curfman McInnes } 884932b0c3eSLois Curfman McInnes } 8850752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 886606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 88790ace30eSBarry Smith } 8883a40ed3dSBarry Smith PetscFunctionReturn(0); 889932b0c3eSLois Curfman McInnes } 890932b0c3eSLois Curfman McInnes 8914a2ae208SSatish Balay #undef __FUNCT__ 8924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 893dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 894f1af5d2fSBarry Smith { 895f1af5d2fSBarry Smith Mat A = (Mat) Aa; 896f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 897*6849ba73SBarry Smith PetscErrorCode ierr; 898*6849ba73SBarry Smith int m = A->m,n = A->n,color,i,j; 89987828ca2SBarry Smith PetscScalar *v = a->v; 900b0a32e0cSBarry Smith PetscViewer viewer; 901b0a32e0cSBarry Smith PetscDraw popup; 902329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 903f3ef73ceSBarry Smith PetscViewerFormat format; 904f1af5d2fSBarry Smith 905f1af5d2fSBarry Smith PetscFunctionBegin; 906f1af5d2fSBarry Smith 907f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 908b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 909b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 910f1af5d2fSBarry Smith 911f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 912fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 913f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 914b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 915f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 916f1af5d2fSBarry Smith x_l = j; 917f1af5d2fSBarry Smith x_r = x_l + 1.0; 918f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 919f1af5d2fSBarry Smith y_l = m - i - 1.0; 920f1af5d2fSBarry Smith y_r = y_l + 1.0; 921f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 922329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 923b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 924329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 925b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 926f1af5d2fSBarry Smith } else { 927f1af5d2fSBarry Smith continue; 928f1af5d2fSBarry Smith } 929f1af5d2fSBarry Smith #else 930f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 931b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 932f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 933b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 934f1af5d2fSBarry Smith } else { 935f1af5d2fSBarry Smith continue; 936f1af5d2fSBarry Smith } 937f1af5d2fSBarry Smith #endif 938b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 939f1af5d2fSBarry Smith } 940f1af5d2fSBarry Smith } 941f1af5d2fSBarry Smith } else { 942f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 943f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 944f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 945f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 946f1af5d2fSBarry Smith } 947b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 948b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 949b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 950f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 951f1af5d2fSBarry Smith x_l = j; 952f1af5d2fSBarry Smith x_r = x_l + 1.0; 953f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 954f1af5d2fSBarry Smith y_l = m - i - 1.0; 955f1af5d2fSBarry Smith y_r = y_l + 1.0; 956b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 957b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 958f1af5d2fSBarry Smith } 959f1af5d2fSBarry Smith } 960f1af5d2fSBarry Smith } 961f1af5d2fSBarry Smith PetscFunctionReturn(0); 962f1af5d2fSBarry Smith } 963f1af5d2fSBarry Smith 9644a2ae208SSatish Balay #undef __FUNCT__ 9654a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 966dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 967f1af5d2fSBarry Smith { 968b0a32e0cSBarry Smith PetscDraw draw; 969f1af5d2fSBarry Smith PetscTruth isnull; 970329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 971dfbe8321SBarry Smith PetscErrorCode ierr; 972f1af5d2fSBarry Smith 973f1af5d2fSBarry Smith PetscFunctionBegin; 974b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 975b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 976f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 977f1af5d2fSBarry Smith 978f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 979273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 980f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 981b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 982b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 983f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 984f1af5d2fSBarry Smith PetscFunctionReturn(0); 985f1af5d2fSBarry Smith } 986f1af5d2fSBarry Smith 9874a2ae208SSatish Balay #undef __FUNCT__ 9884a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 989dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 990932b0c3eSLois Curfman McInnes { 991932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 992dfbe8321SBarry Smith PetscErrorCode ierr; 99332077d6dSBarry Smith PetscTruth issocket,iascii,isbinary,isdraw; 994932b0c3eSLois Curfman McInnes 9953a40ed3dSBarry Smith PetscFunctionBegin; 996b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 99732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 998fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 999fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 10000f5bd95cSBarry Smith 10010f5bd95cSBarry Smith if (issocket) { 1002634064b4SBarry Smith if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA"); 1003b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 100432077d6dSBarry Smith } else if (iascii) { 10053a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10060f5bd95cSBarry Smith } else if (isbinary) { 10073a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1008f1af5d2fSBarry Smith } else if (isdraw) { 1009f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10105cd90555SBarry Smith } else { 1011958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1012932b0c3eSLois Curfman McInnes } 10133a40ed3dSBarry Smith PetscFunctionReturn(0); 1014932b0c3eSLois Curfman McInnes } 1015289bc588SBarry Smith 10164a2ae208SSatish Balay #undef __FUNCT__ 10174a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1018dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1019289bc588SBarry Smith { 1020ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1021dfbe8321SBarry Smith PetscErrorCode ierr; 102290f02eecSBarry Smith 10233a40ed3dSBarry Smith PetscFunctionBegin; 1024aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1025b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 1026a5a9c739SBarry Smith #endif 1027606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 1028606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1029606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 10303a40ed3dSBarry Smith PetscFunctionReturn(0); 1031289bc588SBarry Smith } 1032289bc588SBarry Smith 10334a2ae208SSatish Balay #undef __FUNCT__ 10344a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1035dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout) 1036289bc588SBarry Smith { 1037c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1038*6849ba73SBarry Smith PetscErrorCode ierr; 1039*6849ba73SBarry Smith int k,j,m,n,M; 104087828ca2SBarry Smith PetscScalar *v,tmp; 104148b35521SBarry Smith 10423a40ed3dSBarry Smith PetscFunctionBegin; 10431b807ce4Svictorle v = mat->v; m = A->m; M = mat->lda; n = A->n; 10447c922b88SBarry Smith if (!matout) { /* in place transpose */ 1045a5ce6ee0Svictorle if (m != n) { 1046634064b4SBarry Smith SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1047d64ed03dSBarry Smith } else { 1048d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1049289bc588SBarry Smith for (k=0; k<j; k++) { 10501b807ce4Svictorle tmp = v[j + k*M]; 10511b807ce4Svictorle v[j + k*M] = v[k + j*M]; 10521b807ce4Svictorle v[k + j*M] = tmp; 1053289bc588SBarry Smith } 1054289bc588SBarry Smith } 1055d64ed03dSBarry Smith } 10563a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1057d3e5ee88SLois Curfman McInnes Mat tmat; 1058ec8511deSBarry Smith Mat_SeqDense *tmatd; 105987828ca2SBarry Smith PetscScalar *v2; 1060ea709b57SSatish Balay 10615c5985e7SKris Buschelman ierr = MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);CHKERRQ(ierr); 10625c5985e7SKris Buschelman ierr = MatSetType(tmat,A->type_name);CHKERRQ(ierr); 10635c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1064ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10650de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1066d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10671b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1068d3e5ee88SLois Curfman McInnes } 10696d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10706d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1071d3e5ee88SLois Curfman McInnes *matout = tmat; 107248b35521SBarry Smith } 10733a40ed3dSBarry Smith PetscFunctionReturn(0); 1074289bc588SBarry Smith } 1075289bc588SBarry Smith 10764a2ae208SSatish Balay #undef __FUNCT__ 10774a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1078dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1079289bc588SBarry Smith { 1080c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1081c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1082a43ee2ecSKris Buschelman int i,j; 108387828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 10849ea5d5aeSSatish Balay 10853a40ed3dSBarry Smith PetscFunctionBegin; 1086273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1087273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10881b807ce4Svictorle for (i=0; i<A1->m; i++) { 10891b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 10901b807ce4Svictorle for (j=0; j<A1->n; j++) { 10913a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10921b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 10931b807ce4Svictorle } 1094289bc588SBarry Smith } 109577c4ece6SBarry Smith *flg = PETSC_TRUE; 10963a40ed3dSBarry Smith PetscFunctionReturn(0); 1097289bc588SBarry Smith } 1098289bc588SBarry Smith 10994a2ae208SSatish Balay #undef __FUNCT__ 11004a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1101dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1102289bc588SBarry Smith { 1103c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1104dfbe8321SBarry Smith PetscErrorCode ierr; 1105dfbe8321SBarry Smith int i,n,len; 110687828ca2SBarry Smith PetscScalar *x,zero = 0.0; 110744cd7ae7SLois Curfman McInnes 11083a40ed3dSBarry Smith PetscFunctionBegin; 11097a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 11107a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 11111ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1112273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1113273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 111444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 11151b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1116289bc588SBarry Smith } 11171ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 11183a40ed3dSBarry Smith PetscFunctionReturn(0); 1119289bc588SBarry Smith } 1120289bc588SBarry Smith 11214a2ae208SSatish Balay #undef __FUNCT__ 11224a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1123dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1124289bc588SBarry Smith { 1125c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 112687828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1127dfbe8321SBarry Smith PetscErrorCode ierr; 1128dfbe8321SBarry Smith int i,j,m = A->m,n = A->n; 112955659b69SBarry Smith 11303a40ed3dSBarry Smith PetscFunctionBegin; 113128988994SBarry Smith if (ll) { 11327a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 11331ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1134273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1135da3a660dSBarry Smith for (i=0; i<m; i++) { 1136da3a660dSBarry Smith x = l[i]; 1137da3a660dSBarry Smith v = mat->v + i; 1138da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1139da3a660dSBarry Smith } 11401ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1141b0a32e0cSBarry Smith PetscLogFlops(n*m); 1142da3a660dSBarry Smith } 114328988994SBarry Smith if (rr) { 11447a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 11451ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1146273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1147da3a660dSBarry Smith for (i=0; i<n; i++) { 1148da3a660dSBarry Smith x = r[i]; 1149da3a660dSBarry Smith v = mat->v + i*m; 1150da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1151da3a660dSBarry Smith } 11521ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1153b0a32e0cSBarry Smith PetscLogFlops(n*m); 1154da3a660dSBarry Smith } 11553a40ed3dSBarry Smith PetscFunctionReturn(0); 1156289bc588SBarry Smith } 1157289bc588SBarry Smith 11584a2ae208SSatish Balay #undef __FUNCT__ 11594a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1160dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1161289bc588SBarry Smith { 1162c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 116387828ca2SBarry Smith PetscScalar *v = mat->v; 1164329f5518SBarry Smith PetscReal sum = 0.0; 1165a5ce6ee0Svictorle int lda=mat->lda,m=A->m,i,j; 116655659b69SBarry Smith 11673a40ed3dSBarry Smith PetscFunctionBegin; 1168289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1169a5ce6ee0Svictorle if (lda>m) { 1170a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1171a5ce6ee0Svictorle v = mat->v+j*lda; 1172a5ce6ee0Svictorle for (i=0; i<m; i++) { 1173a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1174a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1175a5ce6ee0Svictorle #else 1176a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1177a5ce6ee0Svictorle #endif 1178a5ce6ee0Svictorle } 1179a5ce6ee0Svictorle } 1180a5ce6ee0Svictorle } else { 1181273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1182aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1183329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1184289bc588SBarry Smith #else 1185289bc588SBarry Smith sum += (*v)*(*v); v++; 1186289bc588SBarry Smith #endif 1187289bc588SBarry Smith } 1188a5ce6ee0Svictorle } 1189064f8208SBarry Smith *nrm = sqrt(sum); 1190b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11913a40ed3dSBarry Smith } else if (type == NORM_1) { 1192064f8208SBarry Smith *nrm = 0.0; 1193273d9f13SBarry Smith for (j=0; j<A->n; j++) { 11941b807ce4Svictorle v = mat->v + j*mat->lda; 1195289bc588SBarry Smith sum = 0.0; 1196273d9f13SBarry Smith for (i=0; i<A->m; i++) { 119733a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1198289bc588SBarry Smith } 1199064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1200289bc588SBarry Smith } 1201b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 12023a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1203064f8208SBarry Smith *nrm = 0.0; 1204273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1205289bc588SBarry Smith v = mat->v + j; 1206289bc588SBarry Smith sum = 0.0; 1207273d9f13SBarry Smith for (i=0; i<A->n; i++) { 12081b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1209289bc588SBarry Smith } 1210064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1211289bc588SBarry Smith } 1212b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 12133a40ed3dSBarry Smith } else { 121429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1215289bc588SBarry Smith } 12163a40ed3dSBarry Smith PetscFunctionReturn(0); 1217289bc588SBarry Smith } 1218289bc588SBarry Smith 12194a2ae208SSatish Balay #undef __FUNCT__ 12204a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1221dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op) 1222289bc588SBarry Smith { 1223c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 122467e560aaSBarry Smith 12253a40ed3dSBarry Smith PetscFunctionBegin; 1226b5a2b587SKris Buschelman switch (op) { 1227b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1228b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1229b5a2b587SKris Buschelman break; 1230b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1231b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1232b5a2b587SKris Buschelman break; 1233b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1234b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1235b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1236b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1237b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1238b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1239b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1240b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1241b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1242b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1243b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1244b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1245b5a2b587SKris Buschelman break; 124677e54ba9SKris Buschelman case MAT_SYMMETRIC: 124777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 12489a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 12499a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 12509a4540c5SBarry Smith case MAT_HERMITIAN: 12519a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 12529a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 12539a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 125477e54ba9SKris Buschelman break; 1255b5a2b587SKris Buschelman default: 125629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 12573a40ed3dSBarry Smith } 12583a40ed3dSBarry Smith PetscFunctionReturn(0); 1259289bc588SBarry Smith } 1260289bc588SBarry Smith 12614a2ae208SSatish Balay #undef __FUNCT__ 12624a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1263dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 12646f0a148fSBarry Smith { 1265ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1266*6849ba73SBarry Smith PetscErrorCode ierr; 1267*6849ba73SBarry Smith int lda=l->lda,m=A->m,j; 12683a40ed3dSBarry Smith 12693a40ed3dSBarry Smith PetscFunctionBegin; 1270a5ce6ee0Svictorle if (lda>m) { 1271a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1272a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1273a5ce6ee0Svictorle } 1274a5ce6ee0Svictorle } else { 127587828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1276a5ce6ee0Svictorle } 12773a40ed3dSBarry Smith PetscFunctionReturn(0); 12786f0a148fSBarry Smith } 12796f0a148fSBarry Smith 12804a2ae208SSatish Balay #undef __FUNCT__ 12814a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 1282dfbe8321SBarry Smith PetscErrorCode MatGetBlockSize_SeqDense(Mat A,int *bs) 12834e220ebcSLois Curfman McInnes { 12843a40ed3dSBarry Smith PetscFunctionBegin; 12854e220ebcSLois Curfman McInnes *bs = 1; 12863a40ed3dSBarry Smith PetscFunctionReturn(0); 12874e220ebcSLois Curfman McInnes } 12884e220ebcSLois Curfman McInnes 12894a2ae208SSatish Balay #undef __FUNCT__ 12904a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1291dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag) 12926f0a148fSBarry Smith { 1293ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1294*6849ba73SBarry Smith PetscErrorCode ierr; 1295*6849ba73SBarry Smith int n = A->n,i,j,N,*rows; 129687828ca2SBarry Smith PetscScalar *slot; 129755659b69SBarry Smith 12983a40ed3dSBarry Smith PetscFunctionBegin; 1299e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 130078b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 13016f0a148fSBarry Smith for (i=0; i<N; i++) { 13026f0a148fSBarry Smith slot = l->v + rows[i]; 13036f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13046f0a148fSBarry Smith } 13056f0a148fSBarry Smith if (diag) { 13066f0a148fSBarry Smith for (i=0; i<N; i++) { 13076f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1308260925b8SBarry Smith *slot = *diag; 13096f0a148fSBarry Smith } 13106f0a148fSBarry Smith } 1311260925b8SBarry Smith ISRestoreIndices(is,&rows); 13123a40ed3dSBarry Smith PetscFunctionReturn(0); 13136f0a148fSBarry Smith } 1314557bce09SLois Curfman McInnes 13154a2ae208SSatish Balay #undef __FUNCT__ 13164a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1317dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 131864e87e97SBarry Smith { 1319c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13203a40ed3dSBarry Smith 13213a40ed3dSBarry Smith PetscFunctionBegin; 132264e87e97SBarry Smith *array = mat->v; 13233a40ed3dSBarry Smith PetscFunctionReturn(0); 132464e87e97SBarry Smith } 13250754003eSLois Curfman McInnes 13264a2ae208SSatish Balay #undef __FUNCT__ 13274a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1328dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1329ff14e315SSatish Balay { 13303a40ed3dSBarry Smith PetscFunctionBegin; 133109b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13323a40ed3dSBarry Smith PetscFunctionReturn(0); 1333ff14e315SSatish Balay } 13340754003eSLois Curfman McInnes 13354a2ae208SSatish Balay #undef __FUNCT__ 13364a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 1337*6849ba73SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 13380754003eSLois Curfman McInnes { 1339c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1340*6849ba73SBarry Smith PetscErrorCode ierr; 1341*6849ba73SBarry Smith int i,j,m = A->m,*irow,*icol,nrows,ncols; 134287828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13430754003eSLois Curfman McInnes Mat newmat; 13440754003eSLois Curfman McInnes 13453a40ed3dSBarry Smith PetscFunctionBegin; 134678b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 134778b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1348e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1349e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13500754003eSLois Curfman McInnes 1351182d2002SSatish Balay /* Check submatrixcall */ 1352182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1353182d2002SSatish Balay int n_cols,n_rows; 1354182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 135529bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1356182d2002SSatish Balay newmat = *B; 1357182d2002SSatish Balay } else { 13580754003eSLois Curfman McInnes /* Create and fill new matrix */ 13595c5985e7SKris Buschelman ierr = MatCreate(A->comm,nrows,ncols,nrows,ncols,&newmat);CHKERRQ(ierr); 13605c5985e7SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 13615c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1362182d2002SSatish Balay } 1363182d2002SSatish Balay 1364182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1365182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1366182d2002SSatish Balay 1367182d2002SSatish Balay for (i=0; i<ncols; i++) { 1368182d2002SSatish Balay av = v + m*icol[i]; 1369182d2002SSatish Balay for (j=0; j<nrows; j++) { 1370182d2002SSatish Balay *bv++ = av[irow[j]]; 13710754003eSLois Curfman McInnes } 13720754003eSLois Curfman McInnes } 1373182d2002SSatish Balay 1374182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13756d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13766d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13770754003eSLois Curfman McInnes 13780754003eSLois Curfman McInnes /* Free work space */ 137978b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 138078b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1381182d2002SSatish Balay *B = newmat; 13823a40ed3dSBarry Smith PetscFunctionReturn(0); 13830754003eSLois Curfman McInnes } 13840754003eSLois Curfman McInnes 13854a2ae208SSatish Balay #undef __FUNCT__ 13864a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 1387dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1388905e6a2fSBarry Smith { 1389*6849ba73SBarry Smith PetscErrorCode ierr; 1390*6849ba73SBarry Smith int i; 1391905e6a2fSBarry Smith 13923a40ed3dSBarry Smith PetscFunctionBegin; 1393905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1394b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1395905e6a2fSBarry Smith } 1396905e6a2fSBarry Smith 1397905e6a2fSBarry Smith for (i=0; i<n; i++) { 13986a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1399905e6a2fSBarry Smith } 14003a40ed3dSBarry Smith PetscFunctionReturn(0); 1401905e6a2fSBarry Smith } 1402905e6a2fSBarry Smith 14034a2ae208SSatish Balay #undef __FUNCT__ 14044a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1405dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14064b0e389bSBarry Smith { 14074b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1408*6849ba73SBarry Smith PetscErrorCode ierr; 1409*6849ba73SBarry Smith int lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j; 14103a40ed3dSBarry Smith 14113a40ed3dSBarry Smith PetscFunctionBegin; 141233f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 141333f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1414cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 14153a40ed3dSBarry Smith PetscFunctionReturn(0); 14163a40ed3dSBarry Smith } 14170dbb7854Svictorle if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1418a5ce6ee0Svictorle if (lda1>m || lda2>m) { 14190dbb7854Svictorle for (j=0; j<n; j++) { 1420a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1421a5ce6ee0Svictorle } 1422a5ce6ee0Svictorle } else { 142387828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1424a5ce6ee0Svictorle } 1425273d9f13SBarry Smith PetscFunctionReturn(0); 1426273d9f13SBarry Smith } 1427273d9f13SBarry Smith 14284a2ae208SSatish Balay #undef __FUNCT__ 14294a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1430dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1431273d9f13SBarry Smith { 1432dfbe8321SBarry Smith PetscErrorCode ierr; 1433273d9f13SBarry Smith 1434273d9f13SBarry Smith PetscFunctionBegin; 1435273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 14363a40ed3dSBarry Smith PetscFunctionReturn(0); 14374b0e389bSBarry Smith } 14384b0e389bSBarry Smith 1439289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1440a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1441905e6a2fSBarry Smith MatGetRow_SeqDense, 1442905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1443905e6a2fSBarry Smith MatMult_SeqDense, 144497304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 14457c922b88SBarry Smith MatMultTranspose_SeqDense, 14467c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1447905e6a2fSBarry Smith MatSolve_SeqDense, 1448905e6a2fSBarry Smith MatSolveAdd_SeqDense, 14497c922b88SBarry Smith MatSolveTranspose_SeqDense, 145097304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense, 1451905e6a2fSBarry Smith MatLUFactor_SeqDense, 1452905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1453ec8511deSBarry Smith MatRelax_SeqDense, 1454ec8511deSBarry Smith MatTranspose_SeqDense, 145597304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1456905e6a2fSBarry Smith MatEqual_SeqDense, 1457905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1458905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1459905e6a2fSBarry Smith MatNorm_SeqDense, 146097304618SKris Buschelman /*20*/ 0, 1461905e6a2fSBarry Smith 0, 1462905e6a2fSBarry Smith 0, 1463905e6a2fSBarry Smith MatSetOption_SeqDense, 1464905e6a2fSBarry Smith MatZeroEntries_SeqDense, 146597304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1466905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1467905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1468905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1469905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 147097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1471273d9f13SBarry Smith 0, 1472905e6a2fSBarry Smith 0, 1473905e6a2fSBarry Smith MatGetArray_SeqDense, 1474905e6a2fSBarry Smith MatRestoreArray_SeqDense, 147597304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1476a5ae1ecdSBarry Smith 0, 1477a5ae1ecdSBarry Smith 0, 1478a5ae1ecdSBarry Smith 0, 1479a5ae1ecdSBarry Smith 0, 148097304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1481a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1482a5ae1ecdSBarry Smith 0, 14834b0e389bSBarry Smith MatGetValues_SeqDense, 1484a5ae1ecdSBarry Smith MatCopy_SeqDense, 148597304618SKris Buschelman /*45*/ 0, 1486a5ae1ecdSBarry Smith MatScale_SeqDense, 1487a5ae1ecdSBarry Smith 0, 1488a5ae1ecdSBarry Smith 0, 1489a5ae1ecdSBarry Smith 0, 149097304618SKris Buschelman /*50*/ MatGetBlockSize_SeqDense, 1491a5ae1ecdSBarry Smith 0, 1492a5ae1ecdSBarry Smith 0, 1493a5ae1ecdSBarry Smith 0, 1494a5ae1ecdSBarry Smith 0, 149597304618SKris Buschelman /*55*/ 0, 1496a5ae1ecdSBarry Smith 0, 1497a5ae1ecdSBarry Smith 0, 1498a5ae1ecdSBarry Smith 0, 1499a5ae1ecdSBarry Smith 0, 150097304618SKris Buschelman /*60*/ 0, 1501e03a110bSBarry Smith MatDestroy_SeqDense, 1502e03a110bSBarry Smith MatView_SeqDense, 150397304618SKris Buschelman MatGetPetscMaps_Petsc, 150497304618SKris Buschelman 0, 150597304618SKris Buschelman /*65*/ 0, 150697304618SKris Buschelman 0, 150797304618SKris Buschelman 0, 150897304618SKris Buschelman 0, 150997304618SKris Buschelman 0, 151097304618SKris Buschelman /*70*/ 0, 151197304618SKris Buschelman 0, 151297304618SKris Buschelman 0, 151397304618SKris Buschelman 0, 151497304618SKris Buschelman 0, 151597304618SKris Buschelman /*75*/ 0, 151697304618SKris Buschelman 0, 151797304618SKris Buschelman 0, 151897304618SKris Buschelman 0, 151997304618SKris Buschelman 0, 152097304618SKris Buschelman /*80*/ 0, 152197304618SKris Buschelman 0, 152297304618SKris Buschelman 0, 152397304618SKris Buschelman 0, 152497304618SKris Buschelman /*85*/ MatLoad_SeqDense}; 152590ace30eSBarry Smith 15264a2ae208SSatish Balay #undef __FUNCT__ 15274a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 15284b828684SBarry Smith /*@C 1529fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1530d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1531d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1532289bc588SBarry Smith 1533db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1534db81eaa0SLois Curfman McInnes 153520563c6bSBarry Smith Input Parameters: 1536db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 15370c775827SLois Curfman McInnes . m - number of rows 153818f449edSLois Curfman McInnes . n - number of columns 1539db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1540dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 154120563c6bSBarry Smith 154220563c6bSBarry Smith Output Parameter: 154344cd7ae7SLois Curfman McInnes . A - the matrix 154420563c6bSBarry Smith 1545b259b22eSLois Curfman McInnes Notes: 154618f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 154718f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1548b4fd4287SBarry Smith set data=PETSC_NULL. 154918f449edSLois Curfman McInnes 1550027ccd11SLois Curfman McInnes Level: intermediate 1551027ccd11SLois Curfman McInnes 1552dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1553d65003e9SLois Curfman McInnes 1554db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 155520563c6bSBarry Smith @*/ 1556dfbe8321SBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1557289bc588SBarry Smith { 1558dfbe8321SBarry Smith PetscErrorCode ierr; 15593b2fbd54SBarry Smith 15603a40ed3dSBarry Smith PetscFunctionBegin; 1561273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1562273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1563273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1564273d9f13SBarry Smith PetscFunctionReturn(0); 1565273d9f13SBarry Smith } 1566273d9f13SBarry Smith 15674a2ae208SSatish Balay #undef __FUNCT__ 15684a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1569273d9f13SBarry Smith /*@C 1570273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1571273d9f13SBarry Smith 1572273d9f13SBarry Smith Collective on MPI_Comm 1573273d9f13SBarry Smith 1574273d9f13SBarry Smith Input Parameters: 1575273d9f13SBarry Smith + A - the matrix 1576273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1577273d9f13SBarry Smith 1578273d9f13SBarry Smith Notes: 1579273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1580273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1581273d9f13SBarry Smith set data=PETSC_NULL. 1582273d9f13SBarry Smith 1583273d9f13SBarry Smith Level: intermediate 1584273d9f13SBarry Smith 1585273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1586273d9f13SBarry Smith 1587273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1588273d9f13SBarry Smith @*/ 1589dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1590273d9f13SBarry Smith { 1591dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1592a23d5eceSKris Buschelman 1593a23d5eceSKris Buschelman PetscFunctionBegin; 1594a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1595a23d5eceSKris Buschelman if (f) { 1596a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1597a23d5eceSKris Buschelman } 1598a23d5eceSKris Buschelman PetscFunctionReturn(0); 1599a23d5eceSKris Buschelman } 1600a23d5eceSKris Buschelman 1601a23d5eceSKris Buschelman EXTERN_C_BEGIN 1602a23d5eceSKris Buschelman #undef __FUNCT__ 1603a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1604dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1605a23d5eceSKris Buschelman { 1606273d9f13SBarry Smith Mat_SeqDense *b; 1607dfbe8321SBarry Smith PetscErrorCode ierr; 1608273d9f13SBarry Smith 1609273d9f13SBarry Smith PetscFunctionBegin; 1610273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1611273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1612273d9f13SBarry Smith if (!data) { 161387828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 161487828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1615273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 161687828ca2SBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1617273d9f13SBarry Smith } else { /* user-allocated storage */ 1618273d9f13SBarry Smith b->v = data; 1619273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1620273d9f13SBarry Smith } 1621273d9f13SBarry Smith PetscFunctionReturn(0); 1622273d9f13SBarry Smith } 1623a23d5eceSKris Buschelman EXTERN_C_END 1624273d9f13SBarry Smith 16251b807ce4Svictorle #undef __FUNCT__ 16261b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 16271b807ce4Svictorle /*@C 16281b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 16291b807ce4Svictorle 16301b807ce4Svictorle Input parameter: 16311b807ce4Svictorle + A - the matrix 16321b807ce4Svictorle - lda - the leading dimension 16331b807ce4Svictorle 16341b807ce4Svictorle Notes: 16351b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 16361b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 16371b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 16381b807ce4Svictorle 16391b807ce4Svictorle Level: intermediate 16401b807ce4Svictorle 16411b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 16421b807ce4Svictorle 16431b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 16441b807ce4Svictorle @*/ 1645dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,int lda) 16461b807ce4Svictorle { 16471b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 16481b807ce4Svictorle PetscFunctionBegin; 1649634064b4SBarry Smith if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %d must be at least matrix dimension %d",lda,B->m); 16501b807ce4Svictorle b->lda = lda; 16511b807ce4Svictorle PetscFunctionReturn(0); 16521b807ce4Svictorle } 16531b807ce4Svictorle 16540bad9183SKris Buschelman /*MC 1655fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 16560bad9183SKris Buschelman 16570bad9183SKris Buschelman Options Database Keys: 16580bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 16590bad9183SKris Buschelman 16600bad9183SKris Buschelman Level: beginner 16610bad9183SKris Buschelman 16620bad9183SKris Buschelman .seealso: MatCreateSeqDense 16630bad9183SKris Buschelman M*/ 16640bad9183SKris Buschelman 1665273d9f13SBarry Smith EXTERN_C_BEGIN 16664a2ae208SSatish Balay #undef __FUNCT__ 16674a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1668dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 1669273d9f13SBarry Smith { 1670273d9f13SBarry Smith Mat_SeqDense *b; 1671dfbe8321SBarry Smith PetscErrorCode ierr; 1672dfbe8321SBarry Smith int size; 1673273d9f13SBarry Smith 1674273d9f13SBarry Smith PetscFunctionBegin; 1675273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 167629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 167755659b69SBarry Smith 1678273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1679273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1680273d9f13SBarry Smith 1681b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1682549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1683549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 168444cd7ae7SLois Curfman McInnes B->factor = 0; 168590f02eecSBarry Smith B->mapping = 0; 1686b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 168744cd7ae7SLois Curfman McInnes B->data = (void*)b; 168818f449edSLois Curfman McInnes 16898a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 16908a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1691a5ae1ecdSBarry Smith 169244cd7ae7SLois Curfman McInnes b->pivots = 0; 1693273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1694273d9f13SBarry Smith b->v = 0; 16951b807ce4Svictorle b->lda = B->m; 16964e220ebcSLois Curfman McInnes 1697a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1698a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 1699a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 17003a40ed3dSBarry Smith PetscFunctionReturn(0); 1701289bc588SBarry Smith } 1702273d9f13SBarry Smith EXTERN_C_END 1703