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; 14a5ce6ee0Svictorle int N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1; 153a40ed3dSBarry Smith 163a40ed3dSBarry Smith PetscFunctionBegin; 17a5ce6ee0Svictorle if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 18a5ce6ee0Svictorle if (ldax>m || lday>m) { 19a5ce6ee0Svictorle for (j=0; j<X->n; j++) { 20a5ce6ee0Svictorle BLaxpy_(&m,alpha,x->v+j*ldax,&one,y->v+j*lday,&one); 21a5ce6ee0Svictorle } 22a5ce6ee0Svictorle } else { 231987afe7SBarry Smith BLaxpy_(&N,alpha,x->v,&one,y->v,&one); 24a5ce6ee0Svictorle } 25b0a32e0cSBarry Smith PetscLogFlops(2*N-1); 263a40ed3dSBarry Smith PetscFunctionReturn(0); 271987afe7SBarry Smith } 281987afe7SBarry Smith 294a2ae208SSatish Balay #undef __FUNCT__ 304a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 318f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 32289bc588SBarry Smith { 334e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 34273d9f13SBarry Smith int i,N = A->m*A->n,count = 0; 3587828ca2SBarry Smith PetscScalar *v = a->v; 363a40ed3dSBarry Smith 373a40ed3dSBarry Smith PetscFunctionBegin; 38289bc588SBarry Smith for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 394e220ebcSLois Curfman McInnes 40273d9f13SBarry Smith info->rows_global = (double)A->m; 41273d9f13SBarry Smith info->columns_global = (double)A->n; 42273d9f13SBarry Smith info->rows_local = (double)A->m; 43273d9f13SBarry Smith info->columns_local = (double)A->n; 444e220ebcSLois Curfman McInnes info->block_size = 1.0; 454e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 464e220ebcSLois Curfman McInnes info->nz_used = (double)count; 474e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 484e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 494e220ebcSLois Curfman McInnes info->mallocs = 0; 504e220ebcSLois Curfman McInnes info->memory = A->mem; 514e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 524e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 534e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 544e220ebcSLois Curfman McInnes 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56289bc588SBarry Smith } 57289bc588SBarry Smith 584a2ae208SSatish Balay #undef __FUNCT__ 594a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 6087828ca2SBarry Smith int MatScale_SeqDense(PetscScalar *alpha,Mat A) 6180cd9d93SLois Curfman McInnes { 62273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 63a5ce6ee0Svictorle int one = 1,lda = a->lda,j,nz; 6480cd9d93SLois Curfman McInnes 653a40ed3dSBarry Smith PetscFunctionBegin; 66a5ce6ee0Svictorle if (lda>A->m) { 67a5ce6ee0Svictorle nz = A->m; 68a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 69a5ce6ee0Svictorle BLscal_(&nz,alpha,a->v+j*lda,&one); 70a5ce6ee0Svictorle } 71a5ce6ee0Svictorle } else { 72273d9f13SBarry Smith nz = A->m*A->n; 7380cd9d93SLois Curfman McInnes BLscal_(&nz,alpha,a->v,&one); 74a5ce6ee0Svictorle } 75b0a32e0cSBarry Smith PetscLogFlops(nz); 763a40ed3dSBarry Smith PetscFunctionReturn(0); 7780cd9d93SLois Curfman McInnes } 7880cd9d93SLois Curfman McInnes 79289bc588SBarry Smith /* ---------------------------------------------------------------*/ 806831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 816831982aSBarry Smith rather than put it in the Mat->row slot.*/ 824a2ae208SSatish Balay #undef __FUNCT__ 834a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 84b380c88cSHong Zhang int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 85289bc588SBarry Smith { 86c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 87b0a32e0cSBarry Smith int info,ierr; 8855659b69SBarry Smith 893a40ed3dSBarry Smith PetscFunctionBegin; 90289bc588SBarry Smith if (!mat->pivots) { 91b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr); 92b0a32e0cSBarry Smith PetscLogObjectMemory(A,A->m*sizeof(int)); 93289bc588SBarry Smith } 94f1af5d2fSBarry Smith A->factor = FACTOR_LU; 95273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 96ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRF) 97ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavilable."); 98ae7cfcebSSatish Balay #else 991b807ce4Svictorle LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info); 10029bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 10129bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 1022ffef20aSBarry Smith #endif 103b0a32e0cSBarry Smith PetscLogFlops((2*A->n*A->n*A->n)/3); 1043a40ed3dSBarry Smith PetscFunctionReturn(0); 105289bc588SBarry Smith } 1066ee01492SSatish Balay 1074a2ae208SSatish Balay #undef __FUNCT__ 1084a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 1095609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11002cad45dSBarry Smith { 11102cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 112a5ce6ee0Svictorle int lda = mat->lda,j,m,ierr; 11302cad45dSBarry Smith Mat newi; 11402cad45dSBarry Smith 1153a40ed3dSBarry Smith PetscFunctionBegin; 116273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr); 1175609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 118a5ce6ee0Svictorle l = (Mat_SeqDense*)newi->data; 119a5ce6ee0Svictorle if (lda>A->m) { 120a5ce6ee0Svictorle m = A->m; 121a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 122a5ce6ee0Svictorle ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 123a5ce6ee0Svictorle } 124a5ce6ee0Svictorle } else { 12587828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 12602cad45dSBarry Smith } 127a5ce6ee0Svictorle } 12802cad45dSBarry Smith newi->assembled = PETSC_TRUE; 12902cad45dSBarry Smith *newmat = newi; 1303a40ed3dSBarry Smith PetscFunctionReturn(0); 13102cad45dSBarry Smith } 13202cad45dSBarry Smith 1334a2ae208SSatish Balay #undef __FUNCT__ 1344a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 135b380c88cSHong Zhang int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 136289bc588SBarry Smith { 1373a40ed3dSBarry Smith int ierr; 1383a40ed3dSBarry Smith 1393a40ed3dSBarry Smith PetscFunctionBegin; 1405609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1413a40ed3dSBarry Smith PetscFunctionReturn(0); 142289bc588SBarry Smith } 1436ee01492SSatish Balay 1444a2ae208SSatish Balay #undef __FUNCT__ 1454a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1468f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 147289bc588SBarry Smith { 14802cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 1491b807ce4Svictorle int lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr; 1503a40ed3dSBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 15202cad45dSBarry Smith /* copy the numerical values */ 1531b807ce4Svictorle if (lda1>m || lda2>m ) { 1541b807ce4Svictorle for (j=0; j<n; j++) { 1551b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar)); CHKERRQ(ierr); 1561b807ce4Svictorle } 1571b807ce4Svictorle } else { 15887828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1591b807ce4Svictorle } 16002cad45dSBarry Smith (*fact)->factor = 0; 161e03a110bSBarry Smith ierr = MatLUFactor(*fact,0,0,PETSC_NULL);CHKERRQ(ierr); 1623a40ed3dSBarry Smith PetscFunctionReturn(0); 163289bc588SBarry Smith } 1646ee01492SSatish Balay 1654a2ae208SSatish Balay #undef __FUNCT__ 1664a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 16715e8a5b3SHong Zhang int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 168289bc588SBarry Smith { 1693a40ed3dSBarry Smith int ierr; 1703a40ed3dSBarry Smith 1713a40ed3dSBarry Smith PetscFunctionBegin; 1723a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1733a40ed3dSBarry Smith PetscFunctionReturn(0); 174289bc588SBarry Smith } 1756ee01492SSatish Balay 1764a2ae208SSatish Balay #undef __FUNCT__ 1774a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 17815e8a5b3SHong Zhang int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 179289bc588SBarry Smith { 180c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 181606d414cSSatish Balay int info,ierr; 18255659b69SBarry Smith 1833a40ed3dSBarry Smith PetscFunctionBegin; 184752f5567SLois Curfman McInnes if (mat->pivots) { 185606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 186b0a32e0cSBarry Smith PetscLogObjectMemory(A,-A->m*sizeof(int)); 187752f5567SLois Curfman McInnes mat->pivots = 0; 188752f5567SLois Curfman McInnes } 189273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 190ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRF) 191ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavilable."); 192ae7cfcebSSatish Balay #else 1931b807ce4Svictorle LApotrf_("L",&A->n,mat->v,&mat->lda,&info); 19429bbc08cSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 1952ffef20aSBarry Smith #endif 196c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 197b0a32e0cSBarry Smith PetscLogFlops((A->n*A->n*A->n)/3); 1983a40ed3dSBarry Smith PetscFunctionReturn(0); 199289bc588SBarry Smith } 200289bc588SBarry Smith 2014a2ae208SSatish Balay #undef __FUNCT__ 2020b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 2030b4b3355SBarry Smith int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 2040b4b3355SBarry Smith { 2050b4b3355SBarry Smith int ierr; 20615e8a5b3SHong Zhang MatFactorInfo info; 2070b4b3355SBarry Smith 2080b4b3355SBarry Smith PetscFunctionBegin; 20915e8a5b3SHong Zhang info.fill = 1.0; 21015e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 2110b4b3355SBarry Smith PetscFunctionReturn(0); 2120b4b3355SBarry Smith } 2130b4b3355SBarry Smith 2140b4b3355SBarry Smith #undef __FUNCT__ 2154a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 2168f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 217289bc588SBarry Smith { 218c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 219a57ad546SLois Curfman McInnes int one = 1,info,ierr; 22087828ca2SBarry Smith PetscScalar *x,*y; 22167e560aaSBarry Smith 2223a40ed3dSBarry Smith PetscFunctionBegin; 223273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2247a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2257a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 22687828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 227c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 228ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 229ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 230ae7cfcebSSatish Balay #else 2311b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2322ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 233ae7cfcebSSatish Balay #endif 2347a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 235ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 236ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 237ae7cfcebSSatish Balay #else 2381b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2392ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 240ae7cfcebSSatish Balay #endif 241289bc588SBarry Smith } 24229bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2437a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2447a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 245b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2463a40ed3dSBarry Smith PetscFunctionReturn(0); 247289bc588SBarry Smith } 2486ee01492SSatish Balay 2494a2ae208SSatish Balay #undef __FUNCT__ 2504a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 2517c922b88SBarry Smith int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 252da3a660dSBarry Smith { 253c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2547a97a34bSBarry Smith int ierr,one = 1,info; 25587828ca2SBarry Smith PetscScalar *x,*y; 25667e560aaSBarry Smith 2573a40ed3dSBarry Smith PetscFunctionBegin; 258273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2597a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2607a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 26187828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 262752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 263da3a660dSBarry Smith if (mat->pivots) { 264ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 265ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 266ae7cfcebSSatish Balay #else 2671b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 2682ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 269ae7cfcebSSatish Balay #endif 2707a97a34bSBarry Smith } else { 271ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 272ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 273ae7cfcebSSatish Balay #else 2741b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 2752ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 276ae7cfcebSSatish Balay #endif 277da3a660dSBarry Smith } 2787a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2797a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 280b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n - A->n); 2813a40ed3dSBarry Smith PetscFunctionReturn(0); 282da3a660dSBarry Smith } 2836ee01492SSatish Balay 2844a2ae208SSatish Balay #undef __FUNCT__ 2854a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 2868f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 287da3a660dSBarry Smith { 288c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2897a97a34bSBarry Smith int ierr,one = 1,info; 29087828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 291da3a660dSBarry Smith Vec tmp = 0; 29267e560aaSBarry Smith 2933a40ed3dSBarry Smith PetscFunctionBegin; 2947a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2957a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 296273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 297da3a660dSBarry Smith if (yy == zz) { 29878b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 299b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 30078b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 301da3a660dSBarry Smith } 30287828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 303752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 304da3a660dSBarry Smith if (mat->pivots) { 305ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 306ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 307ae7cfcebSSatish Balay #else 3081b807ce4Svictorle LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 3092ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 310ae7cfcebSSatish Balay #endif 311a8c6a408SBarry Smith } else { 312ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 313ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 314ae7cfcebSSatish Balay #else 3151b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 3162ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 317ae7cfcebSSatish Balay #endif 318da3a660dSBarry Smith } 319a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 320a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 3217a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3227a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 323b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3243a40ed3dSBarry Smith PetscFunctionReturn(0); 325da3a660dSBarry Smith } 32667e560aaSBarry Smith 3274a2ae208SSatish Balay #undef __FUNCT__ 3284a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 3297c922b88SBarry Smith int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 330da3a660dSBarry Smith { 331c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3326abc6512SBarry Smith int one = 1,info,ierr; 33387828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 334da3a660dSBarry Smith Vec tmp; 33567e560aaSBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 337273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3387a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3397a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 340da3a660dSBarry Smith if (yy == zz) { 34178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 342b0a32e0cSBarry Smith PetscLogObjectParent(A,tmp); 34378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 344da3a660dSBarry Smith } 34587828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 346752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 347da3a660dSBarry Smith if (mat->pivots) { 348ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 349ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavilable."); 350ae7cfcebSSatish Balay #else 3511b807ce4Svictorle LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 3522ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 353ae7cfcebSSatish Balay #endif 3543a40ed3dSBarry Smith } else { 355ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 356ae7cfcebSSatish Balay SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavilable."); 357ae7cfcebSSatish Balay #else 3581b807ce4Svictorle LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 3592ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 360ae7cfcebSSatish Balay #endif 361da3a660dSBarry Smith } 36290f02eecSBarry Smith if (tmp) { 36390f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 36490f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3653a40ed3dSBarry Smith } else { 36690f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 36790f02eecSBarry Smith } 3687a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3697a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 370b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->n); 3713a40ed3dSBarry Smith PetscFunctionReturn(0); 372da3a660dSBarry Smith } 373289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3744a2ae208SSatish Balay #undef __FUNCT__ 3754a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 376329f5518SBarry Smith int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 377c14dc6b6SHong Zhang PetscReal shift,int its,int lits,Vec xx) 378289bc588SBarry Smith { 379c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 38087828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 381273d9f13SBarry Smith int ierr,m = A->m,i; 382aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 383bc1b551cSSatish Balay int o = 1; 384bc1b551cSSatish Balay #endif 385289bc588SBarry Smith 3863a40ed3dSBarry Smith PetscFunctionBegin; 387289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3883a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 389bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 390289bc588SBarry Smith } 3917a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3927a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 393b965ef7fSBarry Smith its = its*lits; 39491723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 395289bc588SBarry Smith while (its--) { 396289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 397289bc588SBarry Smith for (i=0; i<m; i++) { 398aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 399f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 400f1747703SBarry Smith not happy about returning a double complex */ 401f1747703SBarry Smith int _i; 40287828ca2SBarry Smith PetscScalar sum = b[i]; 403f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4043f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 405f1747703SBarry Smith } 406f1747703SBarry Smith xt = sum; 407f1747703SBarry Smith #else 408289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 409f1747703SBarry Smith #endif 41055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 411289bc588SBarry Smith } 412289bc588SBarry Smith } 413289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 414289bc588SBarry Smith for (i=m-1; i>=0; i--) { 415aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 416f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 417f1747703SBarry Smith not happy about returning a double complex */ 418f1747703SBarry Smith int _i; 41987828ca2SBarry Smith PetscScalar sum = b[i]; 420f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4213f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 422f1747703SBarry Smith } 423f1747703SBarry Smith xt = sum; 424f1747703SBarry Smith #else 425289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 426f1747703SBarry Smith #endif 42755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 428289bc588SBarry Smith } 429289bc588SBarry Smith } 430289bc588SBarry Smith } 431600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4327a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4333a40ed3dSBarry Smith PetscFunctionReturn(0); 434289bc588SBarry Smith } 435289bc588SBarry Smith 436289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4374a2ae208SSatish Balay #undef __FUNCT__ 4384a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 4397c922b88SBarry Smith int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 440289bc588SBarry Smith { 441c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 44287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 443ea709b57SSatish Balay int ierr,_One=1; 444ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4453a40ed3dSBarry Smith 4463a40ed3dSBarry Smith PetscFunctionBegin; 447273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4487a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4497a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4501b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4517a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4527a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 453b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->n); 4543a40ed3dSBarry Smith PetscFunctionReturn(0); 455289bc588SBarry Smith } 4566ee01492SSatish Balay 4574a2ae208SSatish Balay #undef __FUNCT__ 4584a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 45944cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 460289bc588SBarry Smith { 461c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 46287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 463329f5518SBarry Smith int ierr,_One=1; 4643a40ed3dSBarry Smith 4653a40ed3dSBarry Smith PetscFunctionBegin; 466273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4677a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4687a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4691b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4707a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4717a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 472b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n - A->m); 4733a40ed3dSBarry Smith PetscFunctionReturn(0); 474289bc588SBarry Smith } 4756ee01492SSatish Balay 4764a2ae208SSatish Balay #undef __FUNCT__ 4774a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 47844cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 479289bc588SBarry Smith { 480c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 48187828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 482329f5518SBarry Smith int ierr,_One=1; 4833a40ed3dSBarry Smith 4843a40ed3dSBarry Smith PetscFunctionBegin; 485273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 486600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4877a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4887a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 4891b807ce4Svictorle LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 4907a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4917a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 492b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 4933a40ed3dSBarry Smith PetscFunctionReturn(0); 494289bc588SBarry Smith } 4956ee01492SSatish Balay 4964a2ae208SSatish Balay #undef __FUNCT__ 4974a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 4987c922b88SBarry Smith int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 499289bc588SBarry Smith { 500c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 50187828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 5027a97a34bSBarry Smith int ierr,_One=1; 50387828ca2SBarry Smith PetscScalar _DOne=1.0; 5043a40ed3dSBarry Smith 5053a40ed3dSBarry Smith PetscFunctionBegin; 506273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 507600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5087a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5097a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5101b807ce4Svictorle LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5117a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5127a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 513b0a32e0cSBarry Smith PetscLogFlops(2*A->m*A->n); 5143a40ed3dSBarry Smith PetscFunctionReturn(0); 515289bc588SBarry Smith } 516289bc588SBarry Smith 517289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5184a2ae208SSatish Balay #undef __FUNCT__ 5194a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 52087828ca2SBarry Smith int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 521289bc588SBarry Smith { 522c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 52387828ca2SBarry Smith PetscScalar *v; 524b0a32e0cSBarry Smith int i,ierr; 52567e560aaSBarry Smith 5263a40ed3dSBarry Smith PetscFunctionBegin; 527273d9f13SBarry Smith *ncols = A->n; 528289bc588SBarry Smith if (cols) { 529b0a32e0cSBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 530273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 531289bc588SBarry Smith } 532289bc588SBarry Smith if (vals) { 53387828ca2SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 534289bc588SBarry Smith v = mat->v + row; 5351b807ce4Svictorle for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;} 536289bc588SBarry Smith } 5373a40ed3dSBarry Smith PetscFunctionReturn(0); 538289bc588SBarry Smith } 5396ee01492SSatish Balay 5404a2ae208SSatish Balay #undef __FUNCT__ 5414a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 54287828ca2SBarry Smith int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 543289bc588SBarry Smith { 544606d414cSSatish Balay int ierr; 545606d414cSSatish Balay PetscFunctionBegin; 546606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 547606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5483a40ed3dSBarry Smith PetscFunctionReturn(0); 549289bc588SBarry Smith } 550289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5514a2ae208SSatish Balay #undef __FUNCT__ 5524a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 5538f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 55487828ca2SBarry Smith int *indexn,PetscScalar *v,InsertMode addv) 555289bc588SBarry Smith { 556c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 557289bc588SBarry Smith int i,j; 558d6dfbf8fSBarry Smith 5593a40ed3dSBarry Smith PetscFunctionBegin; 560289bc588SBarry Smith if (!mat->roworiented) { 561dbb450caSBarry Smith if (addv == INSERT_VALUES) { 562289bc588SBarry Smith for (j=0; j<n; j++) { 5635ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 564aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 56529bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56658804f6dSBarry Smith #endif 567289bc588SBarry Smith for (i=0; i<m; i++) { 5685ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 569aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57029bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 57158804f6dSBarry Smith #endif 5721b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] = *v++; 573289bc588SBarry Smith } 574289bc588SBarry Smith } 5753a40ed3dSBarry Smith } else { 576289bc588SBarry Smith for (j=0; j<n; j++) { 5775ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 578aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 57929bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 58058804f6dSBarry Smith #endif 581289bc588SBarry Smith for (i=0; i<m; i++) { 5825ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 583aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 58429bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 58558804f6dSBarry Smith #endif 5861b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] += *v++; 587289bc588SBarry Smith } 588289bc588SBarry Smith } 589289bc588SBarry Smith } 5903a40ed3dSBarry Smith } else { 591dbb450caSBarry Smith if (addv == INSERT_VALUES) { 592e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 5935ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 594aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 59529bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 59658804f6dSBarry Smith #endif 597e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 5985ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 599aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 60029bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 60158804f6dSBarry Smith #endif 6021b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] = *v++; 603e8d4e0b9SBarry Smith } 604e8d4e0b9SBarry Smith } 6053a40ed3dSBarry Smith } else { 606289bc588SBarry Smith for (i=0; i<m; i++) { 6075ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 608aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 60929bbc08cSBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 61058804f6dSBarry Smith #endif 611289bc588SBarry Smith for (j=0; j<n; j++) { 6125ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 613aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 61429bbc08cSBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 61558804f6dSBarry Smith #endif 6161b807ce4Svictorle mat->v[indexn[j]*mat->lda + indexm[i]] += *v++; 617289bc588SBarry Smith } 618289bc588SBarry Smith } 619289bc588SBarry Smith } 620e8d4e0b9SBarry Smith } 6213a40ed3dSBarry Smith PetscFunctionReturn(0); 622289bc588SBarry Smith } 623e8d4e0b9SBarry Smith 6244a2ae208SSatish Balay #undef __FUNCT__ 6254a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 62687828ca2SBarry Smith int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,PetscScalar *v) 627ae80bb75SLois Curfman McInnes { 628ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 629ae80bb75SLois Curfman McInnes int i,j; 63087828ca2SBarry Smith PetscScalar *vpt = v; 631ae80bb75SLois Curfman McInnes 6323a40ed3dSBarry Smith PetscFunctionBegin; 633ae80bb75SLois Curfman McInnes /* row-oriented output */ 634ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 635ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 6361b807ce4Svictorle *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 637ae80bb75SLois Curfman McInnes } 638ae80bb75SLois Curfman McInnes } 6393a40ed3dSBarry Smith PetscFunctionReturn(0); 640ae80bb75SLois Curfman McInnes } 641ae80bb75SLois Curfman McInnes 642289bc588SBarry Smith /* -----------------------------------------------------------------*/ 643289bc588SBarry Smith 644e090d566SSatish Balay #include "petscsys.h" 64588e32edaSLois Curfman McInnes 646273d9f13SBarry Smith EXTERN_C_BEGIN 6474a2ae208SSatish Balay #undef __FUNCT__ 6484a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 649b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A) 65088e32edaSLois Curfman McInnes { 65188e32edaSLois Curfman McInnes Mat_SeqDense *a; 65288e32edaSLois Curfman McInnes Mat B; 65388e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 65488e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 65587828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 65619bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 65788e32edaSLois Curfman McInnes 6583a40ed3dSBarry Smith PetscFunctionBegin; 659d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 66029bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 661b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6620752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 663552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 66488e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 66588e32edaSLois Curfman McInnes 66690ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 66790ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 66890ace30eSBarry Smith B = *A; 66990ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6704a41a4d3SSatish Balay v = a->v; 6714a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6724a41a4d3SSatish Balay from row major to column major */ 67387828ca2SBarry Smith ierr = PetscMalloc((M*N+1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 67490ace30eSBarry Smith /* read in nonzero values */ 6754a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6764a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6774a41a4d3SSatish Balay for (j=0; j<N; j++) { 6784a41a4d3SSatish Balay for (i=0; i<M; i++) { 6794a41a4d3SSatish Balay *v++ =w[i*N+j]; 6804a41a4d3SSatish Balay } 6814a41a4d3SSatish Balay } 682606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6836d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6846d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68590ace30eSBarry Smith } else { 68688e32edaSLois Curfman McInnes /* read row lengths */ 687b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 6880752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 68988e32edaSLois Curfman McInnes 69088e32edaSLois Curfman McInnes /* create our matrix */ 691b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 69288e32edaSLois Curfman McInnes B = *A; 69388e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 69488e32edaSLois Curfman McInnes v = a->v; 69588e32edaSLois Curfman McInnes 69688e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 697b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 698b0a32e0cSBarry Smith cols = scols; 6990752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 70087828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 701b0a32e0cSBarry Smith vals = svals; 7020752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 70388e32edaSLois Curfman McInnes 70488e32edaSLois Curfman McInnes /* insert into matrix */ 70588e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 70688e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 70788e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 70888e32edaSLois Curfman McInnes } 709606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 710606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 711606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 71288e32edaSLois Curfman McInnes 7136d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7146d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71590ace30eSBarry Smith } 7163a40ed3dSBarry Smith PetscFunctionReturn(0); 71788e32edaSLois Curfman McInnes } 718273d9f13SBarry Smith EXTERN_C_END 71988e32edaSLois Curfman McInnes 720e090d566SSatish Balay #include "petscsys.h" 721932b0c3eSLois Curfman McInnes 7224a2ae208SSatish Balay #undef __FUNCT__ 7234a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 724b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 725289bc588SBarry Smith { 726932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 727fb9695e5SSatish Balay int ierr,i,j; 728fb9695e5SSatish Balay char *name; 72987828ca2SBarry Smith PetscScalar *v; 730f3ef73ceSBarry Smith PetscViewerFormat format; 731932b0c3eSLois Curfman McInnes 7323a40ed3dSBarry Smith PetscFunctionBegin; 733435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 734b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 735456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7363a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 737fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 738b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 739273d9f13SBarry Smith for (i=0; i<A->m; i++) { 74044cd7ae7SLois Curfman McInnes v = a->v + i; 741b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 742273d9f13SBarry Smith for (j=0; j<A->n; j++) { 743aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 744329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 74562b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 746329f5518SBarry Smith } else if (PetscRealPart(*v)) { 74762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7486831982aSBarry Smith } 74980cd9d93SLois Curfman McInnes #else 7506831982aSBarry Smith if (*v) { 75162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 7526831982aSBarry Smith } 75380cd9d93SLois Curfman McInnes #endif 7541b807ce4Svictorle v += a->lda; 75580cd9d93SLois Curfman McInnes } 756b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 75780cd9d93SLois Curfman McInnes } 758b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7593a40ed3dSBarry Smith } else { 760b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 761aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 762ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 76347989497SBarry Smith /* determine if matrix has all real values */ 76447989497SBarry Smith v = a->v; 7651b807ce4Svictorle for (i=0; i<A->m; i++) { 7661b807ce4Svictorle v = a->v + i; 7671b807ce4Svictorle for (j=0; j<A->n; j++) { 768ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 7691b807ce4Svictorle v += a->lda; 7701b807ce4Svictorle } 77147989497SBarry Smith } 77247989497SBarry Smith #endif 773fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7743a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 775b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 776fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 777fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 778ffac6cdbSBarry Smith } 779ffac6cdbSBarry Smith 780273d9f13SBarry Smith for (i=0; i<A->m; i++) { 781932b0c3eSLois Curfman McInnes v = a->v + i; 782273d9f13SBarry Smith for (j=0; j<A->n; j++) { 783aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 78447989497SBarry Smith if (allreal) { 7851b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 78647989497SBarry Smith } else { 7871b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 78847989497SBarry Smith } 789289bc588SBarry Smith #else 7901b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 791289bc588SBarry Smith #endif 7921b807ce4Svictorle v += a->lda; 793289bc588SBarry Smith } 794b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 795289bc588SBarry Smith } 796fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 797b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 798ffac6cdbSBarry Smith } 799b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 800da3a660dSBarry Smith } 801b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8023a40ed3dSBarry Smith PetscFunctionReturn(0); 803289bc588SBarry Smith } 804289bc588SBarry Smith 8054a2ae208SSatish Balay #undef __FUNCT__ 8064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 807b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 808932b0c3eSLois Curfman McInnes { 809932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 810273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 81187828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 812f3ef73ceSBarry Smith PetscViewerFormat format; 813932b0c3eSLois Curfman McInnes 8143a40ed3dSBarry Smith PetscFunctionBegin; 815b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 81690ace30eSBarry Smith 817b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 818fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 81990ace30eSBarry Smith /* store the matrix as a dense matrix */ 820b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 821552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 82290ace30eSBarry Smith col_lens[1] = m; 82390ace30eSBarry Smith col_lens[2] = n; 82490ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8250752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 826606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 82790ace30eSBarry Smith 82890ace30eSBarry Smith /* write out matrix, by rows */ 82987828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 83090ace30eSBarry Smith v = a->v; 83190ace30eSBarry Smith for (i=0; i<m; i++) { 83290ace30eSBarry Smith for (j=0; j<n; j++) { 83390ace30eSBarry Smith vals[i + j*m] = *v++; 83490ace30eSBarry Smith } 83590ace30eSBarry Smith } 8360752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 837606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 83890ace30eSBarry Smith } else { 839b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 840552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 841932b0c3eSLois Curfman McInnes col_lens[1] = m; 842932b0c3eSLois Curfman McInnes col_lens[2] = n; 843932b0c3eSLois Curfman McInnes col_lens[3] = nz; 844932b0c3eSLois Curfman McInnes 845932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 846932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8470752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 848932b0c3eSLois Curfman McInnes 849932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 850932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 851932b0c3eSLois Curfman McInnes ict = 0; 852932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 853932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 854932b0c3eSLois Curfman McInnes } 8550752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 856606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 857932b0c3eSLois Curfman McInnes 858932b0c3eSLois Curfman McInnes /* store nonzero values */ 85987828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 860932b0c3eSLois Curfman McInnes ict = 0; 861932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 862932b0c3eSLois Curfman McInnes v = a->v + i; 863932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 8641b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 865932b0c3eSLois Curfman McInnes } 866932b0c3eSLois Curfman McInnes } 8670752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 868606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 86990ace30eSBarry Smith } 8703a40ed3dSBarry Smith PetscFunctionReturn(0); 871932b0c3eSLois Curfman McInnes } 872932b0c3eSLois Curfman McInnes 8734a2ae208SSatish Balay #undef __FUNCT__ 8744a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 875b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 876f1af5d2fSBarry Smith { 877f1af5d2fSBarry Smith Mat A = (Mat) Aa; 878f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 879fb9695e5SSatish Balay int m = A->m,n = A->n,color,i,j,ierr; 88087828ca2SBarry Smith PetscScalar *v = a->v; 881b0a32e0cSBarry Smith PetscViewer viewer; 882b0a32e0cSBarry Smith PetscDraw popup; 883329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 884f3ef73ceSBarry Smith PetscViewerFormat format; 885f1af5d2fSBarry Smith 886f1af5d2fSBarry Smith PetscFunctionBegin; 887f1af5d2fSBarry Smith 888f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 889b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 890b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 891f1af5d2fSBarry Smith 892f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 893fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 894f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 895b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 896f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 897f1af5d2fSBarry Smith x_l = j; 898f1af5d2fSBarry Smith x_r = x_l + 1.0; 899f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 900f1af5d2fSBarry Smith y_l = m - i - 1.0; 901f1af5d2fSBarry Smith y_r = y_l + 1.0; 902f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 903329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 904b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 905329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 906b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 907f1af5d2fSBarry Smith } else { 908f1af5d2fSBarry Smith continue; 909f1af5d2fSBarry Smith } 910f1af5d2fSBarry Smith #else 911f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 912b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 913f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 914b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 915f1af5d2fSBarry Smith } else { 916f1af5d2fSBarry Smith continue; 917f1af5d2fSBarry Smith } 918f1af5d2fSBarry Smith #endif 919b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 920f1af5d2fSBarry Smith } 921f1af5d2fSBarry Smith } 922f1af5d2fSBarry Smith } else { 923f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 924f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 925f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 926f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 927f1af5d2fSBarry Smith } 928b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 929b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 930b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 931f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 932f1af5d2fSBarry Smith x_l = j; 933f1af5d2fSBarry Smith x_r = x_l + 1.0; 934f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 935f1af5d2fSBarry Smith y_l = m - i - 1.0; 936f1af5d2fSBarry Smith y_r = y_l + 1.0; 937b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 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 } 942f1af5d2fSBarry Smith PetscFunctionReturn(0); 943f1af5d2fSBarry Smith } 944f1af5d2fSBarry Smith 9454a2ae208SSatish Balay #undef __FUNCT__ 9464a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 947b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 948f1af5d2fSBarry Smith { 949b0a32e0cSBarry Smith PetscDraw draw; 950f1af5d2fSBarry Smith PetscTruth isnull; 951329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 952f1af5d2fSBarry Smith int ierr; 953f1af5d2fSBarry Smith 954f1af5d2fSBarry Smith PetscFunctionBegin; 955b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 956b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 957f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 958f1af5d2fSBarry Smith 959f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 960273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 961f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 962b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 963b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 964f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 965f1af5d2fSBarry Smith PetscFunctionReturn(0); 966f1af5d2fSBarry Smith } 967f1af5d2fSBarry Smith 9684a2ae208SSatish Balay #undef __FUNCT__ 9694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 970b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer) 971932b0c3eSLois Curfman McInnes { 972932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 973bcd2baecSBarry Smith int ierr; 974f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 975932b0c3eSLois Curfman McInnes 9763a40ed3dSBarry Smith PetscFunctionBegin; 977b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 978b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 979fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 980fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 9810f5bd95cSBarry Smith 9820f5bd95cSBarry Smith if (issocket) { 9831b807ce4Svictorle if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA"); 984b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9850f5bd95cSBarry Smith } else if (isascii) { 9863a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9870f5bd95cSBarry Smith } else if (isbinary) { 9883a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 989f1af5d2fSBarry Smith } else if (isdraw) { 990f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9915cd90555SBarry Smith } else { 99229bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 993932b0c3eSLois Curfman McInnes } 9943a40ed3dSBarry Smith PetscFunctionReturn(0); 995932b0c3eSLois Curfman McInnes } 996289bc588SBarry Smith 9974a2ae208SSatish Balay #undef __FUNCT__ 9984a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 999c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 1000289bc588SBarry Smith { 1001ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 100290f02eecSBarry Smith int ierr; 100390f02eecSBarry Smith 10043a40ed3dSBarry Smith PetscFunctionBegin; 1005aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1006b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 1007a5a9c739SBarry Smith #endif 1008606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 1009606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1010606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 10113a40ed3dSBarry Smith PetscFunctionReturn(0); 1012289bc588SBarry Smith } 1013289bc588SBarry Smith 10144a2ae208SSatish Balay #undef __FUNCT__ 10154a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 10168f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 1017289bc588SBarry Smith { 1018c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10191b807ce4Svictorle int k,j,m,n,M,ierr; 102087828ca2SBarry Smith PetscScalar *v,tmp; 102148b35521SBarry Smith 10223a40ed3dSBarry Smith PetscFunctionBegin; 10231b807ce4Svictorle v = mat->v; m = A->m; M = mat->lda; n = A->n; 10247c922b88SBarry Smith if (!matout) { /* in place transpose */ 1025a5ce6ee0Svictorle if (m != n) { 1026a5ce6ee0Svictorle SETERRQ(1,"Can not transpose non-square matrix in place"); 1027d64ed03dSBarry Smith } else { 1028d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1029289bc588SBarry Smith for (k=0; k<j; k++) { 10301b807ce4Svictorle tmp = v[j + k*M]; 10311b807ce4Svictorle v[j + k*M] = v[k + j*M]; 10321b807ce4Svictorle v[k + j*M] = tmp; 1033289bc588SBarry Smith } 1034289bc588SBarry Smith } 1035d64ed03dSBarry Smith } 10363a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1037d3e5ee88SLois Curfman McInnes Mat tmat; 1038ec8511deSBarry Smith Mat_SeqDense *tmatd; 103987828ca2SBarry Smith PetscScalar *v2; 1040ea709b57SSatish Balay 1041273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1042ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10430de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1044d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10451b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1046d3e5ee88SLois Curfman McInnes } 10476d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10486d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1049d3e5ee88SLois Curfman McInnes *matout = tmat; 105048b35521SBarry Smith } 10513a40ed3dSBarry Smith PetscFunctionReturn(0); 1052289bc588SBarry Smith } 1053289bc588SBarry Smith 10544a2ae208SSatish Balay #undef __FUNCT__ 10554a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 10568f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1057289bc588SBarry Smith { 1058c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1059c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 10601b807ce4Svictorle int ierr,i,j; 106187828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1062273d9f13SBarry Smith PetscTruth flag; 10639ea5d5aeSSatish Balay 10643a40ed3dSBarry Smith PetscFunctionBegin; 1065273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)A2,MATSEQDENSE,&flag);CHKERRQ(ierr); 1066273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_SUP,"Matrices should be of type MATSEQDENSE"); 1067273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1068273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10691b807ce4Svictorle for (i=0; i<A1->m; i++) { 10701b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 10711b807ce4Svictorle for (j=0; j<A1->n; j++) { 10723a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10731b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 10741b807ce4Svictorle } 1075289bc588SBarry Smith } 107677c4ece6SBarry Smith *flg = PETSC_TRUE; 10773a40ed3dSBarry Smith PetscFunctionReturn(0); 1078289bc588SBarry Smith } 1079289bc588SBarry Smith 10804a2ae208SSatish Balay #undef __FUNCT__ 10814a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 10828f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1083289bc588SBarry Smith { 1084c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10857a97a34bSBarry Smith int ierr,i,n,len; 108687828ca2SBarry Smith PetscScalar *x,zero = 0.0; 108744cd7ae7SLois Curfman McInnes 10883a40ed3dSBarry Smith PetscFunctionBegin; 10897a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10907a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1091600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1092273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1093273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 109444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 10951b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1096289bc588SBarry Smith } 10977a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10983a40ed3dSBarry Smith PetscFunctionReturn(0); 1099289bc588SBarry Smith } 1100289bc588SBarry Smith 11014a2ae208SSatish Balay #undef __FUNCT__ 11024a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 11038f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1104289bc588SBarry Smith { 1105c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 110687828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1107273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 110855659b69SBarry Smith 11093a40ed3dSBarry Smith PetscFunctionBegin; 111028988994SBarry Smith if (ll) { 11117a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1112600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1113273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1114da3a660dSBarry Smith for (i=0; i<m; i++) { 1115da3a660dSBarry Smith x = l[i]; 1116da3a660dSBarry Smith v = mat->v + i; 1117da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1118da3a660dSBarry Smith } 11197a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1120b0a32e0cSBarry Smith PetscLogFlops(n*m); 1121da3a660dSBarry Smith } 112228988994SBarry Smith if (rr) { 11237a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1124600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1125273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1126da3a660dSBarry Smith for (i=0; i<n; i++) { 1127da3a660dSBarry Smith x = r[i]; 1128da3a660dSBarry Smith v = mat->v + i*m; 1129da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1130da3a660dSBarry Smith } 11317a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1132b0a32e0cSBarry Smith PetscLogFlops(n*m); 1133da3a660dSBarry Smith } 11343a40ed3dSBarry Smith PetscFunctionReturn(0); 1135289bc588SBarry Smith } 1136289bc588SBarry Smith 11374a2ae208SSatish Balay #undef __FUNCT__ 11384a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1139064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1140289bc588SBarry Smith { 1141c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 114287828ca2SBarry Smith PetscScalar *v = mat->v; 1143329f5518SBarry Smith PetscReal sum = 0.0; 1144a5ce6ee0Svictorle int lda=mat->lda,m=A->m,i,j; 114555659b69SBarry Smith 11463a40ed3dSBarry Smith PetscFunctionBegin; 1147289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1148a5ce6ee0Svictorle if (lda>m) { 1149a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1150a5ce6ee0Svictorle v = mat->v+j*lda; 1151a5ce6ee0Svictorle for (i=0; i<m; i++) { 1152a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1153a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1154a5ce6ee0Svictorle #else 1155a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1156a5ce6ee0Svictorle #endif 1157a5ce6ee0Svictorle } 1158a5ce6ee0Svictorle } 1159a5ce6ee0Svictorle } else { 1160273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1161aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1162329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1163289bc588SBarry Smith #else 1164289bc588SBarry Smith sum += (*v)*(*v); v++; 1165289bc588SBarry Smith #endif 1166289bc588SBarry Smith } 1167a5ce6ee0Svictorle } 1168064f8208SBarry Smith *nrm = sqrt(sum); 1169b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11703a40ed3dSBarry Smith } else if (type == NORM_1) { 1171064f8208SBarry Smith *nrm = 0.0; 1172273d9f13SBarry Smith for (j=0; j<A->n; j++) { 11731b807ce4Svictorle v = mat->v + j*mat->lda; 1174289bc588SBarry Smith sum = 0.0; 1175273d9f13SBarry Smith for (i=0; i<A->m; i++) { 117633a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1177289bc588SBarry Smith } 1178064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1179289bc588SBarry Smith } 1180b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11813a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1182064f8208SBarry Smith *nrm = 0.0; 1183273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1184289bc588SBarry Smith v = mat->v + j; 1185289bc588SBarry Smith sum = 0.0; 1186273d9f13SBarry Smith for (i=0; i<A->n; i++) { 11871b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1188289bc588SBarry Smith } 1189064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1190289bc588SBarry Smith } 1191b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11923a40ed3dSBarry Smith } else { 119329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1194289bc588SBarry Smith } 11953a40ed3dSBarry Smith PetscFunctionReturn(0); 1196289bc588SBarry Smith } 1197289bc588SBarry Smith 11984a2ae208SSatish Balay #undef __FUNCT__ 11994a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 12008f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1201289bc588SBarry Smith { 1202c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 120367e560aaSBarry Smith 12043a40ed3dSBarry Smith PetscFunctionBegin; 1205b5a2b587SKris Buschelman switch (op) { 1206b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1207b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1208b5a2b587SKris Buschelman break; 1209b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1210b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1211b5a2b587SKris Buschelman break; 1212b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1213b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1214b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1215b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1216b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1217b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1218b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1219b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1220b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1221b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1222b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1223b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1224b5a2b587SKris Buschelman break; 1225b5a2b587SKris Buschelman default: 122629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 12273a40ed3dSBarry Smith } 12283a40ed3dSBarry Smith PetscFunctionReturn(0); 1229289bc588SBarry Smith } 1230289bc588SBarry Smith 12314a2ae208SSatish Balay #undef __FUNCT__ 12324a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 12338f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 12346f0a148fSBarry Smith { 1235ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1236a5ce6ee0Svictorle int lda=l->lda,m=A->m,j, ierr; 12373a40ed3dSBarry Smith 12383a40ed3dSBarry Smith PetscFunctionBegin; 1239a5ce6ee0Svictorle if (lda>m) { 1240a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1241a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar)); CHKERRQ(ierr); 1242a5ce6ee0Svictorle } 1243a5ce6ee0Svictorle } else { 124487828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1245a5ce6ee0Svictorle } 12463a40ed3dSBarry Smith PetscFunctionReturn(0); 12476f0a148fSBarry Smith } 12486f0a148fSBarry Smith 12494a2ae208SSatish Balay #undef __FUNCT__ 12504a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 12518f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 12524e220ebcSLois Curfman McInnes { 12533a40ed3dSBarry Smith PetscFunctionBegin; 12544e220ebcSLois Curfman McInnes *bs = 1; 12553a40ed3dSBarry Smith PetscFunctionReturn(0); 12564e220ebcSLois Curfman McInnes } 12574e220ebcSLois Curfman McInnes 12584a2ae208SSatish Balay #undef __FUNCT__ 12594a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 126087828ca2SBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,PetscScalar *diag) 12616f0a148fSBarry Smith { 1262ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1263273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 126487828ca2SBarry Smith PetscScalar *slot; 126555659b69SBarry Smith 12663a40ed3dSBarry Smith PetscFunctionBegin; 1267e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 126878b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12696f0a148fSBarry Smith for (i=0; i<N; i++) { 12706f0a148fSBarry Smith slot = l->v + rows[i]; 12716f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 12726f0a148fSBarry Smith } 12736f0a148fSBarry Smith if (diag) { 12746f0a148fSBarry Smith for (i=0; i<N; i++) { 12756f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1276260925b8SBarry Smith *slot = *diag; 12776f0a148fSBarry Smith } 12786f0a148fSBarry Smith } 1279260925b8SBarry Smith ISRestoreIndices(is,&rows); 12803a40ed3dSBarry Smith PetscFunctionReturn(0); 12816f0a148fSBarry Smith } 1282557bce09SLois Curfman McInnes 12834a2ae208SSatish Balay #undef __FUNCT__ 12844a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 128587828ca2SBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar **array) 128664e87e97SBarry Smith { 1287c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12883a40ed3dSBarry Smith 12893a40ed3dSBarry Smith PetscFunctionBegin; 129064e87e97SBarry Smith *array = mat->v; 12913a40ed3dSBarry Smith PetscFunctionReturn(0); 129264e87e97SBarry Smith } 12930754003eSLois Curfman McInnes 12944a2ae208SSatish Balay #undef __FUNCT__ 12954a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 129687828ca2SBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar **array) 1297ff14e315SSatish Balay { 12983a40ed3dSBarry Smith PetscFunctionBegin; 129909b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13003a40ed3dSBarry Smith PetscFunctionReturn(0); 1301ff14e315SSatish Balay } 13020754003eSLois Curfman McInnes 13034a2ae208SSatish Balay #undef __FUNCT__ 13044a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 13057b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 13060754003eSLois Curfman McInnes { 1307c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1308273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 130987828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13100754003eSLois Curfman McInnes Mat newmat; 13110754003eSLois Curfman McInnes 13123a40ed3dSBarry Smith PetscFunctionBegin; 131378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 131478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1315e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1316e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13170754003eSLois Curfman McInnes 1318182d2002SSatish Balay /* Check submatrixcall */ 1319182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1320182d2002SSatish Balay int n_cols,n_rows; 1321182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 132229bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1323182d2002SSatish Balay newmat = *B; 1324182d2002SSatish Balay } else { 13250754003eSLois Curfman McInnes /* Create and fill new matrix */ 1326b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1327182d2002SSatish Balay } 1328182d2002SSatish Balay 1329182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1330182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1331182d2002SSatish Balay 1332182d2002SSatish Balay for (i=0; i<ncols; i++) { 1333182d2002SSatish Balay av = v + m*icol[i]; 1334182d2002SSatish Balay for (j=0; j<nrows; j++) { 1335182d2002SSatish Balay *bv++ = av[irow[j]]; 13360754003eSLois Curfman McInnes } 13370754003eSLois Curfman McInnes } 1338182d2002SSatish Balay 1339182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13406d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13416d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13420754003eSLois Curfman McInnes 13430754003eSLois Curfman McInnes /* Free work space */ 134478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 134578b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1346182d2002SSatish Balay *B = newmat; 13473a40ed3dSBarry Smith PetscFunctionReturn(0); 13480754003eSLois Curfman McInnes } 13490754003eSLois Curfman McInnes 13504a2ae208SSatish Balay #undef __FUNCT__ 13514a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 13527b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 1353905e6a2fSBarry Smith { 1354905e6a2fSBarry Smith int ierr,i; 1355905e6a2fSBarry Smith 13563a40ed3dSBarry Smith PetscFunctionBegin; 1357905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1358b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1359905e6a2fSBarry Smith } 1360905e6a2fSBarry Smith 1361905e6a2fSBarry Smith for (i=0; i<n; i++) { 13626a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1363905e6a2fSBarry Smith } 13643a40ed3dSBarry Smith PetscFunctionReturn(0); 1365905e6a2fSBarry Smith } 1366905e6a2fSBarry Smith 13674a2ae208SSatish Balay #undef __FUNCT__ 13684a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1369cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 13704b0e389bSBarry Smith { 13714b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1372a5ce6ee0Svictorle int lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr; 1373273d9f13SBarry Smith PetscTruth flag; 13743a40ed3dSBarry Smith 13753a40ed3dSBarry Smith PetscFunctionBegin; 1376273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&flag);CHKERRQ(ierr); 1377273d9f13SBarry Smith if (!flag) { 1378cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 13793a40ed3dSBarry Smith PetscFunctionReturn(0); 13803a40ed3dSBarry Smith } 13810dbb7854Svictorle if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1382a5ce6ee0Svictorle if (lda1>m || lda2>m) { 13830dbb7854Svictorle for (j=0; j<n; j++) { 1384a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1385a5ce6ee0Svictorle } 1386a5ce6ee0Svictorle } else { 138787828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1388a5ce6ee0Svictorle } 1389273d9f13SBarry Smith PetscFunctionReturn(0); 1390273d9f13SBarry Smith } 1391273d9f13SBarry Smith 13924a2ae208SSatish Balay #undef __FUNCT__ 13934a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1394273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1395273d9f13SBarry Smith { 1396273d9f13SBarry Smith int ierr; 1397273d9f13SBarry Smith 1398273d9f13SBarry Smith PetscFunctionBegin; 1399273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 14003a40ed3dSBarry Smith PetscFunctionReturn(0); 14014b0e389bSBarry Smith } 14024b0e389bSBarry Smith 1403289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1404a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1405905e6a2fSBarry Smith MatGetRow_SeqDense, 1406905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1407905e6a2fSBarry Smith MatMult_SeqDense, 1408905e6a2fSBarry Smith MatMultAdd_SeqDense, 14097c922b88SBarry Smith MatMultTranspose_SeqDense, 14107c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1411905e6a2fSBarry Smith MatSolve_SeqDense, 1412905e6a2fSBarry Smith MatSolveAdd_SeqDense, 14137c922b88SBarry Smith MatSolveTranspose_SeqDense, 14147c922b88SBarry Smith MatSolveTransposeAdd_SeqDense, 1415905e6a2fSBarry Smith MatLUFactor_SeqDense, 1416905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1417ec8511deSBarry Smith MatRelax_SeqDense, 1418ec8511deSBarry Smith MatTranspose_SeqDense, 1419905e6a2fSBarry Smith MatGetInfo_SeqDense, 1420905e6a2fSBarry Smith MatEqual_SeqDense, 1421905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1422905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1423905e6a2fSBarry Smith MatNorm_SeqDense, 1424905e6a2fSBarry Smith 0, 1425905e6a2fSBarry Smith 0, 1426905e6a2fSBarry Smith 0, 1427905e6a2fSBarry Smith MatSetOption_SeqDense, 1428905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1429905e6a2fSBarry Smith MatZeroRows_SeqDense, 1430905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1431905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1432905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1433905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1434273d9f13SBarry Smith MatSetUpPreallocation_SeqDense, 1435273d9f13SBarry Smith 0, 1436905e6a2fSBarry Smith 0, 1437905e6a2fSBarry Smith MatGetArray_SeqDense, 1438905e6a2fSBarry Smith MatRestoreArray_SeqDense, 14395609ef8eSBarry Smith MatDuplicate_SeqDense, 1440a5ae1ecdSBarry Smith 0, 1441a5ae1ecdSBarry Smith 0, 1442a5ae1ecdSBarry Smith 0, 1443a5ae1ecdSBarry Smith 0, 1444a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1445a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1446a5ae1ecdSBarry Smith 0, 14474b0e389bSBarry Smith MatGetValues_SeqDense, 1448a5ae1ecdSBarry Smith MatCopy_SeqDense, 1449a5ae1ecdSBarry Smith 0, 1450a5ae1ecdSBarry Smith MatScale_SeqDense, 1451a5ae1ecdSBarry Smith 0, 1452a5ae1ecdSBarry Smith 0, 1453a5ae1ecdSBarry Smith 0, 1454a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1455a5ae1ecdSBarry Smith 0, 1456a5ae1ecdSBarry Smith 0, 1457a5ae1ecdSBarry Smith 0, 1458a5ae1ecdSBarry Smith 0, 1459a5ae1ecdSBarry Smith 0, 1460a5ae1ecdSBarry Smith 0, 1461a5ae1ecdSBarry Smith 0, 1462a5ae1ecdSBarry Smith 0, 1463a5ae1ecdSBarry Smith 0, 1464a5ae1ecdSBarry Smith 0, 1465e03a110bSBarry Smith MatDestroy_SeqDense, 1466e03a110bSBarry Smith MatView_SeqDense, 14678a124369SBarry Smith MatGetPetscMaps_Petsc}; 146890ace30eSBarry Smith 14694a2ae208SSatish Balay #undef __FUNCT__ 14704a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 14714b828684SBarry Smith /*@C 1472fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1473d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1474d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1475289bc588SBarry Smith 1476db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1477db81eaa0SLois Curfman McInnes 147820563c6bSBarry Smith Input Parameters: 1479db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 14800c775827SLois Curfman McInnes . m - number of rows 148118f449edSLois Curfman McInnes . n - number of columns 1482db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1483dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 148420563c6bSBarry Smith 148520563c6bSBarry Smith Output Parameter: 148644cd7ae7SLois Curfman McInnes . A - the matrix 148720563c6bSBarry Smith 1488b259b22eSLois Curfman McInnes Notes: 148918f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 149018f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1491b4fd4287SBarry Smith set data=PETSC_NULL. 149218f449edSLois Curfman McInnes 1493027ccd11SLois Curfman McInnes Level: intermediate 1494027ccd11SLois Curfman McInnes 1495dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1496d65003e9SLois Curfman McInnes 1497db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 149820563c6bSBarry Smith @*/ 149987828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1500289bc588SBarry Smith { 1501273d9f13SBarry Smith int ierr; 15023b2fbd54SBarry Smith 15033a40ed3dSBarry Smith PetscFunctionBegin; 1504273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1505273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1506273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1507273d9f13SBarry Smith PetscFunctionReturn(0); 1508273d9f13SBarry Smith } 1509273d9f13SBarry Smith 15104a2ae208SSatish Balay #undef __FUNCT__ 15114a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1512273d9f13SBarry Smith /*@C 1513273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1514273d9f13SBarry Smith 1515273d9f13SBarry Smith Collective on MPI_Comm 1516273d9f13SBarry Smith 1517273d9f13SBarry Smith Input Parameters: 1518273d9f13SBarry Smith + A - the matrix 1519273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1520273d9f13SBarry Smith 1521273d9f13SBarry Smith Notes: 1522273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1523273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1524273d9f13SBarry Smith set data=PETSC_NULL. 1525273d9f13SBarry Smith 1526273d9f13SBarry Smith Level: intermediate 1527273d9f13SBarry Smith 1528273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1529273d9f13SBarry Smith 1530273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1531273d9f13SBarry Smith @*/ 153287828ca2SBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar *data) 1533273d9f13SBarry Smith { 1534*a23d5eceSKris Buschelman int ierr,(*f)(Mat,PetscScalar*); 1535*a23d5eceSKris Buschelman 1536*a23d5eceSKris Buschelman PetscFunctionBegin; 1537*a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1538*a23d5eceSKris Buschelman if (f) { 1539*a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1540*a23d5eceSKris Buschelman } 1541*a23d5eceSKris Buschelman PetscFunctionReturn(0); 1542*a23d5eceSKris Buschelman } 1543*a23d5eceSKris Buschelman 1544*a23d5eceSKris Buschelman EXTERN_C_BEGIN 1545*a23d5eceSKris Buschelman #undef __FUNCT__ 1546*a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1547*a23d5eceSKris Buschelman int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1548*a23d5eceSKris Buschelman { 1549273d9f13SBarry Smith Mat_SeqDense *b; 1550273d9f13SBarry Smith int ierr; 1551273d9f13SBarry Smith PetscTruth flg2; 1552273d9f13SBarry Smith 1553273d9f13SBarry Smith PetscFunctionBegin; 1554273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1555273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1556273d9f13SBarry Smith if (!data) { 155787828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 155887828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1559273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 156087828ca2SBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1561273d9f13SBarry Smith } else { /* user-allocated storage */ 1562273d9f13SBarry Smith b->v = data; 1563273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1564273d9f13SBarry Smith } 1565273d9f13SBarry Smith PetscFunctionReturn(0); 1566273d9f13SBarry Smith } 1567*a23d5eceSKris Buschelman EXTERN_C_END 1568273d9f13SBarry Smith 15691b807ce4Svictorle #undef __FUNCT__ 15701b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 15711b807ce4Svictorle /*@C 15721b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 15731b807ce4Svictorle 15741b807ce4Svictorle Input parameter: 15751b807ce4Svictorle + A - the matrix 15761b807ce4Svictorle - lda - the leading dimension 15771b807ce4Svictorle 15781b807ce4Svictorle Notes: 15791b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 15801b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 15811b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 15821b807ce4Svictorle 15831b807ce4Svictorle Level: intermediate 15841b807ce4Svictorle 15851b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 15861b807ce4Svictorle 15871b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 15881b807ce4Svictorle @*/ 15891b807ce4Svictorle int MatSeqDenseSetLDA(Mat B,int lda) 15901b807ce4Svictorle { 15911b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 15921b807ce4Svictorle PetscFunctionBegin; 15931b807ce4Svictorle if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension"); 15941b807ce4Svictorle b->lda = lda; 15951b807ce4Svictorle PetscFunctionReturn(0); 15961b807ce4Svictorle } 15971b807ce4Svictorle 1598273d9f13SBarry Smith EXTERN_C_BEGIN 15994a2ae208SSatish Balay #undef __FUNCT__ 16004a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1601273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1602273d9f13SBarry Smith { 1603273d9f13SBarry Smith Mat_SeqDense *b; 1604273d9f13SBarry Smith int ierr,size; 1605273d9f13SBarry Smith 1606273d9f13SBarry Smith PetscFunctionBegin; 1607273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 160829bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 160955659b69SBarry Smith 1610273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1611273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1612273d9f13SBarry Smith 1613b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1614549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1615549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 161644cd7ae7SLois Curfman McInnes B->factor = 0; 161790f02eecSBarry Smith B->mapping = 0; 1618b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 161944cd7ae7SLois Curfman McInnes B->data = (void*)b; 162018f449edSLois Curfman McInnes 16218a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 16228a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1623a5ae1ecdSBarry Smith 162444cd7ae7SLois Curfman McInnes b->pivots = 0; 1625273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1626273d9f13SBarry Smith b->v = 0; 16271b807ce4Svictorle b->lda = B->m; 16284e220ebcSLois Curfman McInnes 1629*a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1630*a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 1631*a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 16323a40ed3dSBarry Smith PetscFunctionReturn(0); 1633289bc588SBarry Smith } 1634273d9f13SBarry Smith EXTERN_C_END 1635