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" 11268466fbSBarry Smith int MatAXPY_SeqDense(const 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++) { 20268466fbSBarry Smith BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one); 21a5ce6ee0Svictorle } 22a5ce6ee0Svictorle } else { 23268466fbSBarry Smith BLaxpy_(&N,(PetscScalar*)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" 60268466fbSBarry Smith int MatScale_SeqDense(const 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++) { 69268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one); 70a5ce6ee0Svictorle } 71a5ce6ee0Svictorle } else { 72273d9f13SBarry Smith nz = A->m*A->n; 73268466fbSBarry Smith BLscal_(&nz,(PetscScalar*)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" 553f15d580aSBarry Smith int MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv) 554289bc588SBarry Smith { 555c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 556cddbea37SSatish Balay int i,j,idx=0; 557d6dfbf8fSBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 559289bc588SBarry Smith if (!mat->roworiented) { 560dbb450caSBarry Smith if (addv == INSERT_VALUES) { 561289bc588SBarry Smith for (j=0; j<n; j++) { 562cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 563aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 564590ac198SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 56558804f6dSBarry Smith #endif 566289bc588SBarry Smith for (i=0; i<m; i++) { 567cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 568aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 5695c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 57058804f6dSBarry Smith #endif 571cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 572289bc588SBarry Smith } 573289bc588SBarry Smith } 5743a40ed3dSBarry Smith } else { 575289bc588SBarry Smith for (j=0; j<n; j++) { 576cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 577aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 578590ac198SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 57958804f6dSBarry Smith #endif 580289bc588SBarry Smith for (i=0; i<m; i++) { 581cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 582aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 5835c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 58458804f6dSBarry Smith #endif 585cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 586289bc588SBarry Smith } 587289bc588SBarry Smith } 588289bc588SBarry Smith } 5893a40ed3dSBarry Smith } else { 590dbb450caSBarry Smith if (addv == INSERT_VALUES) { 591e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 592cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 593aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 5945c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 59558804f6dSBarry Smith #endif 596e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 597cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 598aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 5995c8e7597SSatish Balay if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 60058804f6dSBarry Smith #endif 601cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 602e8d4e0b9SBarry Smith } 603e8d4e0b9SBarry Smith } 6043a40ed3dSBarry Smith } else { 605289bc588SBarry Smith for (i=0; i<m; i++) { 606cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 607aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6085c8e7597SSatish Balay if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 60958804f6dSBarry Smith #endif 610289bc588SBarry Smith for (j=0; j<n; j++) { 611cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 612aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 6135c8e7597SSatish Balay if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 61458804f6dSBarry Smith #endif 615cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 616289bc588SBarry Smith } 617289bc588SBarry Smith } 618289bc588SBarry Smith } 619e8d4e0b9SBarry Smith } 6203a40ed3dSBarry Smith PetscFunctionReturn(0); 621289bc588SBarry Smith } 622e8d4e0b9SBarry Smith 6234a2ae208SSatish Balay #undef __FUNCT__ 6244a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 625f15d580aSBarry Smith int MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[]) 626ae80bb75SLois Curfman McInnes { 627ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 628ae80bb75SLois Curfman McInnes int i,j; 62987828ca2SBarry Smith PetscScalar *vpt = v; 630ae80bb75SLois Curfman McInnes 6313a40ed3dSBarry Smith PetscFunctionBegin; 632ae80bb75SLois Curfman McInnes /* row-oriented output */ 633ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 634ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 6351b807ce4Svictorle *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 636ae80bb75SLois Curfman McInnes } 637ae80bb75SLois Curfman McInnes } 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 639ae80bb75SLois Curfman McInnes } 640ae80bb75SLois Curfman McInnes 641289bc588SBarry Smith /* -----------------------------------------------------------------*/ 642289bc588SBarry Smith 643e090d566SSatish Balay #include "petscsys.h" 64488e32edaSLois Curfman McInnes 6454a2ae208SSatish Balay #undef __FUNCT__ 6464a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 647b0a32e0cSBarry Smith int MatLoad_SeqDense(PetscViewer viewer,MatType type,Mat *A) 64888e32edaSLois Curfman McInnes { 64988e32edaSLois Curfman McInnes Mat_SeqDense *a; 65088e32edaSLois Curfman McInnes Mat B; 65188e32edaSLois Curfman McInnes int *scols,i,j,nz,ierr,fd,header[4],size; 65288e32edaSLois Curfman McInnes int *rowlengths = 0,M,N,*cols; 65387828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 65419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 65588e32edaSLois Curfman McInnes 6563a40ed3dSBarry Smith PetscFunctionBegin; 657d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 65829bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 659b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6600752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 661552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 66288e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 66388e32edaSLois Curfman McInnes 66490ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 66590ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 66690ace30eSBarry Smith B = *A; 66790ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6684a41a4d3SSatish Balay v = a->v; 6694a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6704a41a4d3SSatish Balay from row major to column major */ 671b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 67290ace30eSBarry Smith /* read in nonzero values */ 6734a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6744a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6754a41a4d3SSatish Balay for (j=0; j<N; j++) { 6764a41a4d3SSatish Balay for (i=0; i<M; i++) { 6774a41a4d3SSatish Balay *v++ =w[i*N+j]; 6784a41a4d3SSatish Balay } 6794a41a4d3SSatish Balay } 680606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6816d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6826d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68390ace30eSBarry Smith } else { 68488e32edaSLois Curfman McInnes /* read row lengths */ 685b0a32e0cSBarry Smith ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 6860752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 68788e32edaSLois Curfman McInnes 68888e32edaSLois Curfman McInnes /* create our matrix */ 689b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 69088e32edaSLois Curfman McInnes B = *A; 69188e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 69288e32edaSLois Curfman McInnes v = a->v; 69388e32edaSLois Curfman McInnes 69488e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 695b0a32e0cSBarry Smith ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 696b0a32e0cSBarry Smith cols = scols; 6970752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 69887828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 699b0a32e0cSBarry Smith vals = svals; 7000752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 70188e32edaSLois Curfman McInnes 70288e32edaSLois Curfman McInnes /* insert into matrix */ 70388e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 70488e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 70588e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 70688e32edaSLois Curfman McInnes } 707606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 708606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 709606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 71088e32edaSLois Curfman McInnes 7116d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7126d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71390ace30eSBarry Smith } 7143a40ed3dSBarry Smith PetscFunctionReturn(0); 71588e32edaSLois Curfman McInnes } 71688e32edaSLois Curfman McInnes 717e090d566SSatish Balay #include "petscsys.h" 718932b0c3eSLois Curfman McInnes 7194a2ae208SSatish Balay #undef __FUNCT__ 7204a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 721b0a32e0cSBarry Smith static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 722289bc588SBarry Smith { 723932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 724fb9695e5SSatish Balay int ierr,i,j; 725fb9695e5SSatish Balay char *name; 72687828ca2SBarry Smith PetscScalar *v; 727f3ef73ceSBarry Smith PetscViewerFormat format; 728932b0c3eSLois Curfman McInnes 7293a40ed3dSBarry Smith PetscFunctionBegin; 730435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 731b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 732456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7333a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 734fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 735b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 736273d9f13SBarry Smith for (i=0; i<A->m; i++) { 73744cd7ae7SLois Curfman McInnes v = a->v + i; 738b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 739273d9f13SBarry Smith for (j=0; j<A->n; j++) { 740aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 741329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 74262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 743329f5518SBarry Smith } else if (PetscRealPart(*v)) { 74462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7456831982aSBarry Smith } 74680cd9d93SLois Curfman McInnes #else 7476831982aSBarry Smith if (*v) { 74862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 7496831982aSBarry Smith } 75080cd9d93SLois Curfman McInnes #endif 7511b807ce4Svictorle v += a->lda; 75280cd9d93SLois Curfman McInnes } 753b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 75480cd9d93SLois Curfman McInnes } 755b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7563a40ed3dSBarry Smith } else { 757b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 758aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 759ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 76047989497SBarry Smith /* determine if matrix has all real values */ 76147989497SBarry Smith v = a->v; 762201f5f94SSatish Balay for (i=0; i<A->m*A->n; i++) { 763ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 76447989497SBarry Smith } 76547989497SBarry Smith #endif 766fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7673a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 768b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 769fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 770fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 771ffac6cdbSBarry Smith } 772ffac6cdbSBarry Smith 773273d9f13SBarry Smith for (i=0; i<A->m; i++) { 774932b0c3eSLois Curfman McInnes v = a->v + i; 775273d9f13SBarry Smith for (j=0; j<A->n; j++) { 776aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 77747989497SBarry Smith if (allreal) { 7781b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 77947989497SBarry Smith } else { 7801b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 78147989497SBarry Smith } 782289bc588SBarry Smith #else 7831b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 784289bc588SBarry Smith #endif 7851b807ce4Svictorle v += a->lda; 786289bc588SBarry Smith } 787b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 788289bc588SBarry Smith } 789fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 790b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 791ffac6cdbSBarry Smith } 792b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 793da3a660dSBarry Smith } 794b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 7953a40ed3dSBarry Smith PetscFunctionReturn(0); 796289bc588SBarry Smith } 797289bc588SBarry Smith 7984a2ae208SSatish Balay #undef __FUNCT__ 7994a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 800b0a32e0cSBarry Smith static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 801932b0c3eSLois Curfman McInnes { 802932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 803273d9f13SBarry Smith int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 80487828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 805f3ef73ceSBarry Smith PetscViewerFormat format; 806932b0c3eSLois Curfman McInnes 8073a40ed3dSBarry Smith PetscFunctionBegin; 808b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 80990ace30eSBarry Smith 810b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 811fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 81290ace30eSBarry Smith /* store the matrix as a dense matrix */ 813b0a32e0cSBarry Smith ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 814552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 81590ace30eSBarry Smith col_lens[1] = m; 81690ace30eSBarry Smith col_lens[2] = n; 81790ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8180752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 819606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 82090ace30eSBarry Smith 82190ace30eSBarry Smith /* write out matrix, by rows */ 82287828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 82390ace30eSBarry Smith v = a->v; 82490ace30eSBarry Smith for (i=0; i<m; i++) { 82590ace30eSBarry Smith for (j=0; j<n; j++) { 82690ace30eSBarry Smith vals[i + j*m] = *v++; 82790ace30eSBarry Smith } 82890ace30eSBarry Smith } 8290752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 830606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 83190ace30eSBarry Smith } else { 832b0a32e0cSBarry Smith ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 833552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 834932b0c3eSLois Curfman McInnes col_lens[1] = m; 835932b0c3eSLois Curfman McInnes col_lens[2] = n; 836932b0c3eSLois Curfman McInnes col_lens[3] = nz; 837932b0c3eSLois Curfman McInnes 838932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 839932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8400752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 841932b0c3eSLois Curfman McInnes 842932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 843932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 844932b0c3eSLois Curfman McInnes ict = 0; 845932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 846932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 847932b0c3eSLois Curfman McInnes } 8480752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 849606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 850932b0c3eSLois Curfman McInnes 851932b0c3eSLois Curfman McInnes /* store nonzero values */ 85287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 853932b0c3eSLois Curfman McInnes ict = 0; 854932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 855932b0c3eSLois Curfman McInnes v = a->v + i; 856932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 8571b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 858932b0c3eSLois Curfman McInnes } 859932b0c3eSLois Curfman McInnes } 8600752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 861606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 86290ace30eSBarry Smith } 8633a40ed3dSBarry Smith PetscFunctionReturn(0); 864932b0c3eSLois Curfman McInnes } 865932b0c3eSLois Curfman McInnes 8664a2ae208SSatish Balay #undef __FUNCT__ 8674a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 868b0a32e0cSBarry Smith int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 869f1af5d2fSBarry Smith { 870f1af5d2fSBarry Smith Mat A = (Mat) Aa; 871f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 872fb9695e5SSatish Balay int m = A->m,n = A->n,color,i,j,ierr; 87387828ca2SBarry Smith PetscScalar *v = a->v; 874b0a32e0cSBarry Smith PetscViewer viewer; 875b0a32e0cSBarry Smith PetscDraw popup; 876329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 877f3ef73ceSBarry Smith PetscViewerFormat format; 878f1af5d2fSBarry Smith 879f1af5d2fSBarry Smith PetscFunctionBegin; 880f1af5d2fSBarry Smith 881f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 882b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 883b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 884f1af5d2fSBarry Smith 885f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 886fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 887f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 888b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 889f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 890f1af5d2fSBarry Smith x_l = j; 891f1af5d2fSBarry Smith x_r = x_l + 1.0; 892f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 893f1af5d2fSBarry Smith y_l = m - i - 1.0; 894f1af5d2fSBarry Smith y_r = y_l + 1.0; 895f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 896329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 897b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 898329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 899b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 900f1af5d2fSBarry Smith } else { 901f1af5d2fSBarry Smith continue; 902f1af5d2fSBarry Smith } 903f1af5d2fSBarry Smith #else 904f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 905b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 906f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 907b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 908f1af5d2fSBarry Smith } else { 909f1af5d2fSBarry Smith continue; 910f1af5d2fSBarry Smith } 911f1af5d2fSBarry Smith #endif 912b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 913f1af5d2fSBarry Smith } 914f1af5d2fSBarry Smith } 915f1af5d2fSBarry Smith } else { 916f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 917f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 918f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 919f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 920f1af5d2fSBarry Smith } 921b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 922b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 923b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 924f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 925f1af5d2fSBarry Smith x_l = j; 926f1af5d2fSBarry Smith x_r = x_l + 1.0; 927f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 928f1af5d2fSBarry Smith y_l = m - i - 1.0; 929f1af5d2fSBarry Smith y_r = y_l + 1.0; 930b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 931b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 932f1af5d2fSBarry Smith } 933f1af5d2fSBarry Smith } 934f1af5d2fSBarry Smith } 935f1af5d2fSBarry Smith PetscFunctionReturn(0); 936f1af5d2fSBarry Smith } 937f1af5d2fSBarry Smith 9384a2ae208SSatish Balay #undef __FUNCT__ 9394a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 940b0a32e0cSBarry Smith int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 941f1af5d2fSBarry Smith { 942b0a32e0cSBarry Smith PetscDraw draw; 943f1af5d2fSBarry Smith PetscTruth isnull; 944329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 945f1af5d2fSBarry Smith int ierr; 946f1af5d2fSBarry Smith 947f1af5d2fSBarry Smith PetscFunctionBegin; 948b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 949b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 950f1af5d2fSBarry Smith if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 951f1af5d2fSBarry Smith 952f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 953273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 954f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 955b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 956b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 957f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 958f1af5d2fSBarry Smith PetscFunctionReturn(0); 959f1af5d2fSBarry Smith } 960f1af5d2fSBarry Smith 9614a2ae208SSatish Balay #undef __FUNCT__ 9624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 963b0a32e0cSBarry Smith int MatView_SeqDense(Mat A,PetscViewer viewer) 964932b0c3eSLois Curfman McInnes { 965932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 966bcd2baecSBarry Smith int ierr; 967f1af5d2fSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 968932b0c3eSLois Curfman McInnes 9693a40ed3dSBarry Smith PetscFunctionBegin; 970b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 971b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 972fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 973fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 9740f5bd95cSBarry Smith 9750f5bd95cSBarry Smith if (issocket) { 9761b807ce4Svictorle if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA"); 977b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 9780f5bd95cSBarry Smith } else if (isascii) { 9793a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 9800f5bd95cSBarry Smith } else if (isbinary) { 9813a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 982f1af5d2fSBarry Smith } else if (isdraw) { 983f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 9845cd90555SBarry Smith } else { 98529bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 986932b0c3eSLois Curfman McInnes } 9873a40ed3dSBarry Smith PetscFunctionReturn(0); 988932b0c3eSLois Curfman McInnes } 989289bc588SBarry Smith 9904a2ae208SSatish Balay #undef __FUNCT__ 9914a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 992c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 993289bc588SBarry Smith { 994ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 99590f02eecSBarry Smith int ierr; 99690f02eecSBarry Smith 9973a40ed3dSBarry Smith PetscFunctionBegin; 998aa482453SBarry Smith #if defined(PETSC_USE_LOG) 999b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 1000a5a9c739SBarry Smith #endif 1001606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 1002606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1003606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 10043a40ed3dSBarry Smith PetscFunctionReturn(0); 1005289bc588SBarry Smith } 1006289bc588SBarry Smith 10074a2ae208SSatish Balay #undef __FUNCT__ 10084a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 10098f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 1010289bc588SBarry Smith { 1011c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10121b807ce4Svictorle int k,j,m,n,M,ierr; 101387828ca2SBarry Smith PetscScalar *v,tmp; 101448b35521SBarry Smith 10153a40ed3dSBarry Smith PetscFunctionBegin; 10161b807ce4Svictorle v = mat->v; m = A->m; M = mat->lda; n = A->n; 10177c922b88SBarry Smith if (!matout) { /* in place transpose */ 1018a5ce6ee0Svictorle if (m != n) { 1019a5ce6ee0Svictorle SETERRQ(1,"Can not transpose non-square matrix in place"); 1020d64ed03dSBarry Smith } else { 1021d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1022289bc588SBarry Smith for (k=0; k<j; k++) { 10231b807ce4Svictorle tmp = v[j + k*M]; 10241b807ce4Svictorle v[j + k*M] = v[k + j*M]; 10251b807ce4Svictorle v[k + j*M] = tmp; 1026289bc588SBarry Smith } 1027289bc588SBarry Smith } 1028d64ed03dSBarry Smith } 10293a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1030d3e5ee88SLois Curfman McInnes Mat tmat; 1031ec8511deSBarry Smith Mat_SeqDense *tmatd; 103287828ca2SBarry Smith PetscScalar *v2; 1033ea709b57SSatish Balay 1034273d9f13SBarry Smith ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1035ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10360de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1037d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10381b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1039d3e5ee88SLois Curfman McInnes } 10406d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10416d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1042d3e5ee88SLois Curfman McInnes *matout = tmat; 104348b35521SBarry Smith } 10443a40ed3dSBarry Smith PetscFunctionReturn(0); 1045289bc588SBarry Smith } 1046289bc588SBarry Smith 10474a2ae208SSatish Balay #undef __FUNCT__ 10484a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 10498f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1050289bc588SBarry Smith { 1051c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1052c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1053a43ee2ecSKris Buschelman int i,j; 105487828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 10559ea5d5aeSSatish Balay 10563a40ed3dSBarry Smith PetscFunctionBegin; 1057273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1058273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10591b807ce4Svictorle for (i=0; i<A1->m; i++) { 10601b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 10611b807ce4Svictorle for (j=0; j<A1->n; j++) { 10623a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10631b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 10641b807ce4Svictorle } 1065289bc588SBarry Smith } 106677c4ece6SBarry Smith *flg = PETSC_TRUE; 10673a40ed3dSBarry Smith PetscFunctionReturn(0); 1068289bc588SBarry Smith } 1069289bc588SBarry Smith 10704a2ae208SSatish Balay #undef __FUNCT__ 10714a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 10728f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 1073289bc588SBarry Smith { 1074c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10757a97a34bSBarry Smith int ierr,i,n,len; 107687828ca2SBarry Smith PetscScalar *x,zero = 0.0; 107744cd7ae7SLois Curfman McInnes 10783a40ed3dSBarry Smith PetscFunctionBegin; 10797a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 10807a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1081600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1082273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1083273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 108444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 10851b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1086289bc588SBarry Smith } 10877a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 10883a40ed3dSBarry Smith PetscFunctionReturn(0); 1089289bc588SBarry Smith } 1090289bc588SBarry Smith 10914a2ae208SSatish Balay #undef __FUNCT__ 10924a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 10938f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1094289bc588SBarry Smith { 1095c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 109687828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1097273d9f13SBarry Smith int ierr,i,j,m = A->m,n = A->n; 109855659b69SBarry Smith 10993a40ed3dSBarry Smith PetscFunctionBegin; 110028988994SBarry Smith if (ll) { 11017a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1102600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1103273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1104da3a660dSBarry Smith for (i=0; i<m; i++) { 1105da3a660dSBarry Smith x = l[i]; 1106da3a660dSBarry Smith v = mat->v + i; 1107da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1108da3a660dSBarry Smith } 11097a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1110b0a32e0cSBarry Smith PetscLogFlops(n*m); 1111da3a660dSBarry Smith } 111228988994SBarry Smith if (rr) { 11137a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1114600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1115273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1116da3a660dSBarry Smith for (i=0; i<n; i++) { 1117da3a660dSBarry Smith x = r[i]; 1118da3a660dSBarry Smith v = mat->v + i*m; 1119da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1120da3a660dSBarry Smith } 11217a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1122b0a32e0cSBarry Smith PetscLogFlops(n*m); 1123da3a660dSBarry Smith } 11243a40ed3dSBarry Smith PetscFunctionReturn(0); 1125289bc588SBarry Smith } 1126289bc588SBarry Smith 11274a2ae208SSatish Balay #undef __FUNCT__ 11284a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1129064f8208SBarry Smith int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1130289bc588SBarry Smith { 1131c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 113287828ca2SBarry Smith PetscScalar *v = mat->v; 1133329f5518SBarry Smith PetscReal sum = 0.0; 1134a5ce6ee0Svictorle int lda=mat->lda,m=A->m,i,j; 113555659b69SBarry Smith 11363a40ed3dSBarry Smith PetscFunctionBegin; 1137289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1138a5ce6ee0Svictorle if (lda>m) { 1139a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1140a5ce6ee0Svictorle v = mat->v+j*lda; 1141a5ce6ee0Svictorle for (i=0; i<m; i++) { 1142a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1143a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1144a5ce6ee0Svictorle #else 1145a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1146a5ce6ee0Svictorle #endif 1147a5ce6ee0Svictorle } 1148a5ce6ee0Svictorle } 1149a5ce6ee0Svictorle } else { 1150273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1151aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1152329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1153289bc588SBarry Smith #else 1154289bc588SBarry Smith sum += (*v)*(*v); v++; 1155289bc588SBarry Smith #endif 1156289bc588SBarry Smith } 1157a5ce6ee0Svictorle } 1158064f8208SBarry Smith *nrm = sqrt(sum); 1159b0a32e0cSBarry Smith PetscLogFlops(2*A->n*A->m); 11603a40ed3dSBarry Smith } else if (type == NORM_1) { 1161064f8208SBarry Smith *nrm = 0.0; 1162273d9f13SBarry Smith for (j=0; j<A->n; j++) { 11631b807ce4Svictorle v = mat->v + j*mat->lda; 1164289bc588SBarry Smith sum = 0.0; 1165273d9f13SBarry Smith for (i=0; i<A->m; i++) { 116633a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1167289bc588SBarry Smith } 1168064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1169289bc588SBarry Smith } 1170b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11713a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1172064f8208SBarry Smith *nrm = 0.0; 1173273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1174289bc588SBarry Smith v = mat->v + j; 1175289bc588SBarry Smith sum = 0.0; 1176273d9f13SBarry Smith for (i=0; i<A->n; i++) { 11771b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1178289bc588SBarry Smith } 1179064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1180289bc588SBarry Smith } 1181b0a32e0cSBarry Smith PetscLogFlops(A->n*A->m); 11823a40ed3dSBarry Smith } else { 118329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1184289bc588SBarry Smith } 11853a40ed3dSBarry Smith PetscFunctionReturn(0); 1186289bc588SBarry Smith } 1187289bc588SBarry Smith 11884a2ae208SSatish Balay #undef __FUNCT__ 11894a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 11908f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1191289bc588SBarry Smith { 1192c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 119367e560aaSBarry Smith 11943a40ed3dSBarry Smith PetscFunctionBegin; 1195b5a2b587SKris Buschelman switch (op) { 1196b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1197b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1198b5a2b587SKris Buschelman break; 1199b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1200b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1201b5a2b587SKris Buschelman break; 1202b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1203b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1204b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1205b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1206b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1207b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1208b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1209b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1210b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1211b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1212b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1213b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1214b5a2b587SKris Buschelman break; 1215*77e54ba9SKris Buschelman case MAT_SYMMETRIC: 1216*77e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1217*77e54ba9SKris Buschelman break; 1218b5a2b587SKris Buschelman default: 121929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 12203a40ed3dSBarry Smith } 12213a40ed3dSBarry Smith PetscFunctionReturn(0); 1222289bc588SBarry Smith } 1223289bc588SBarry Smith 12244a2ae208SSatish Balay #undef __FUNCT__ 12254a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 12268f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 12276f0a148fSBarry Smith { 1228ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1229a5ce6ee0Svictorle int lda=l->lda,m=A->m,j, ierr; 12303a40ed3dSBarry Smith 12313a40ed3dSBarry Smith PetscFunctionBegin; 1232a5ce6ee0Svictorle if (lda>m) { 1233a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1234a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar)); CHKERRQ(ierr); 1235a5ce6ee0Svictorle } 1236a5ce6ee0Svictorle } else { 123787828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1238a5ce6ee0Svictorle } 12393a40ed3dSBarry Smith PetscFunctionReturn(0); 12406f0a148fSBarry Smith } 12416f0a148fSBarry Smith 12424a2ae208SSatish Balay #undef __FUNCT__ 12434a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqDense" 12448f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 12454e220ebcSLois Curfman McInnes { 12463a40ed3dSBarry Smith PetscFunctionBegin; 12474e220ebcSLois Curfman McInnes *bs = 1; 12483a40ed3dSBarry Smith PetscFunctionReturn(0); 12494e220ebcSLois Curfman McInnes } 12504e220ebcSLois Curfman McInnes 12514a2ae208SSatish Balay #undef __FUNCT__ 12524a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1253268466fbSBarry Smith int MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag) 12546f0a148fSBarry Smith { 1255ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1256273d9f13SBarry Smith int n = A->n,i,j,ierr,N,*rows; 125787828ca2SBarry Smith PetscScalar *slot; 125855659b69SBarry Smith 12593a40ed3dSBarry Smith PetscFunctionBegin; 1260e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 126178b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12626f0a148fSBarry Smith for (i=0; i<N; i++) { 12636f0a148fSBarry Smith slot = l->v + rows[i]; 12646f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 12656f0a148fSBarry Smith } 12666f0a148fSBarry Smith if (diag) { 12676f0a148fSBarry Smith for (i=0; i<N; i++) { 12686f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1269260925b8SBarry Smith *slot = *diag; 12706f0a148fSBarry Smith } 12716f0a148fSBarry Smith } 1272260925b8SBarry Smith ISRestoreIndices(is,&rows); 12733a40ed3dSBarry Smith PetscFunctionReturn(0); 12746f0a148fSBarry Smith } 1275557bce09SLois Curfman McInnes 12764a2ae208SSatish Balay #undef __FUNCT__ 12774a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 12784e7234bfSBarry Smith int MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 127964e87e97SBarry Smith { 1280c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 12813a40ed3dSBarry Smith 12823a40ed3dSBarry Smith PetscFunctionBegin; 128364e87e97SBarry Smith *array = mat->v; 12843a40ed3dSBarry Smith PetscFunctionReturn(0); 128564e87e97SBarry Smith } 12860754003eSLois Curfman McInnes 12874a2ae208SSatish Balay #undef __FUNCT__ 12884a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 12894e7234bfSBarry Smith int MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1290ff14e315SSatish Balay { 12913a40ed3dSBarry Smith PetscFunctionBegin; 129209b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 12933a40ed3dSBarry Smith PetscFunctionReturn(0); 1294ff14e315SSatish Balay } 12950754003eSLois Curfman McInnes 12964a2ae208SSatish Balay #undef __FUNCT__ 12974a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 12987b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 12990754003eSLois Curfman McInnes { 1300c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1301273d9f13SBarry Smith int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 130287828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13030754003eSLois Curfman McInnes Mat newmat; 13040754003eSLois Curfman McInnes 13053a40ed3dSBarry Smith PetscFunctionBegin; 130678b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 130778b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1308e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1309e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13100754003eSLois Curfman McInnes 1311182d2002SSatish Balay /* Check submatrixcall */ 1312182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1313182d2002SSatish Balay int n_cols,n_rows; 1314182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 131529bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1316182d2002SSatish Balay newmat = *B; 1317182d2002SSatish Balay } else { 13180754003eSLois Curfman McInnes /* Create and fill new matrix */ 1319b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1320182d2002SSatish Balay } 1321182d2002SSatish Balay 1322182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1323182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1324182d2002SSatish Balay 1325182d2002SSatish Balay for (i=0; i<ncols; i++) { 1326182d2002SSatish Balay av = v + m*icol[i]; 1327182d2002SSatish Balay for (j=0; j<nrows; j++) { 1328182d2002SSatish Balay *bv++ = av[irow[j]]; 13290754003eSLois Curfman McInnes } 13300754003eSLois Curfman McInnes } 1331182d2002SSatish Balay 1332182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13336d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13346d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13350754003eSLois Curfman McInnes 13360754003eSLois Curfman McInnes /* Free work space */ 133778b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 133878b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1339182d2002SSatish Balay *B = newmat; 13403a40ed3dSBarry Smith PetscFunctionReturn(0); 13410754003eSLois Curfman McInnes } 13420754003eSLois Curfman McInnes 13434a2ae208SSatish Balay #undef __FUNCT__ 13444a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 1345268466fbSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1346905e6a2fSBarry Smith { 1347905e6a2fSBarry Smith int ierr,i; 1348905e6a2fSBarry Smith 13493a40ed3dSBarry Smith PetscFunctionBegin; 1350905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1351b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1352905e6a2fSBarry Smith } 1353905e6a2fSBarry Smith 1354905e6a2fSBarry Smith for (i=0; i<n; i++) { 13556a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1356905e6a2fSBarry Smith } 13573a40ed3dSBarry Smith PetscFunctionReturn(0); 1358905e6a2fSBarry Smith } 1359905e6a2fSBarry Smith 13604a2ae208SSatish Balay #undef __FUNCT__ 13614a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1362cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 13634b0e389bSBarry Smith { 13644b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1365a5ce6ee0Svictorle int lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr; 13663a40ed3dSBarry Smith 13673a40ed3dSBarry Smith PetscFunctionBegin; 136833f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 136933f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1370cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 13713a40ed3dSBarry Smith PetscFunctionReturn(0); 13723a40ed3dSBarry Smith } 13730dbb7854Svictorle if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1374a5ce6ee0Svictorle if (lda1>m || lda2>m) { 13750dbb7854Svictorle for (j=0; j<n; j++) { 1376a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1377a5ce6ee0Svictorle } 1378a5ce6ee0Svictorle } else { 137987828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1380a5ce6ee0Svictorle } 1381273d9f13SBarry Smith PetscFunctionReturn(0); 1382273d9f13SBarry Smith } 1383273d9f13SBarry Smith 13844a2ae208SSatish Balay #undef __FUNCT__ 13854a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1386273d9f13SBarry Smith int MatSetUpPreallocation_SeqDense(Mat A) 1387273d9f13SBarry Smith { 1388273d9f13SBarry Smith int ierr; 1389273d9f13SBarry Smith 1390273d9f13SBarry Smith PetscFunctionBegin; 1391273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 13923a40ed3dSBarry Smith PetscFunctionReturn(0); 13934b0e389bSBarry Smith } 13944b0e389bSBarry Smith 1395289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1396a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1397905e6a2fSBarry Smith MatGetRow_SeqDense, 1398905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1399905e6a2fSBarry Smith MatMult_SeqDense, 140097304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 14017c922b88SBarry Smith MatMultTranspose_SeqDense, 14027c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1403905e6a2fSBarry Smith MatSolve_SeqDense, 1404905e6a2fSBarry Smith MatSolveAdd_SeqDense, 14057c922b88SBarry Smith MatSolveTranspose_SeqDense, 140697304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense, 1407905e6a2fSBarry Smith MatLUFactor_SeqDense, 1408905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1409ec8511deSBarry Smith MatRelax_SeqDense, 1410ec8511deSBarry Smith MatTranspose_SeqDense, 141197304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1412905e6a2fSBarry Smith MatEqual_SeqDense, 1413905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1414905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1415905e6a2fSBarry Smith MatNorm_SeqDense, 141697304618SKris Buschelman /*20*/ 0, 1417905e6a2fSBarry Smith 0, 1418905e6a2fSBarry Smith 0, 1419905e6a2fSBarry Smith MatSetOption_SeqDense, 1420905e6a2fSBarry Smith MatZeroEntries_SeqDense, 142197304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1422905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1423905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1424905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1425905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 142697304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1427273d9f13SBarry Smith 0, 1428905e6a2fSBarry Smith 0, 1429905e6a2fSBarry Smith MatGetArray_SeqDense, 1430905e6a2fSBarry Smith MatRestoreArray_SeqDense, 143197304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1432a5ae1ecdSBarry Smith 0, 1433a5ae1ecdSBarry Smith 0, 1434a5ae1ecdSBarry Smith 0, 1435a5ae1ecdSBarry Smith 0, 143697304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1437a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1438a5ae1ecdSBarry Smith 0, 14394b0e389bSBarry Smith MatGetValues_SeqDense, 1440a5ae1ecdSBarry Smith MatCopy_SeqDense, 144197304618SKris Buschelman /*45*/ 0, 1442a5ae1ecdSBarry Smith MatScale_SeqDense, 1443a5ae1ecdSBarry Smith 0, 1444a5ae1ecdSBarry Smith 0, 1445a5ae1ecdSBarry Smith 0, 144697304618SKris Buschelman /*50*/ MatGetBlockSize_SeqDense, 1447a5ae1ecdSBarry Smith 0, 1448a5ae1ecdSBarry Smith 0, 1449a5ae1ecdSBarry Smith 0, 1450a5ae1ecdSBarry Smith 0, 145197304618SKris Buschelman /*55*/ 0, 1452a5ae1ecdSBarry Smith 0, 1453a5ae1ecdSBarry Smith 0, 1454a5ae1ecdSBarry Smith 0, 1455a5ae1ecdSBarry Smith 0, 145697304618SKris Buschelman /*60*/ 0, 1457e03a110bSBarry Smith MatDestroy_SeqDense, 1458e03a110bSBarry Smith MatView_SeqDense, 145997304618SKris Buschelman MatGetPetscMaps_Petsc, 146097304618SKris Buschelman 0, 146197304618SKris Buschelman /*65*/ 0, 146297304618SKris Buschelman 0, 146397304618SKris Buschelman 0, 146497304618SKris Buschelman 0, 146597304618SKris Buschelman 0, 146697304618SKris Buschelman /*70*/ 0, 146797304618SKris Buschelman 0, 146897304618SKris Buschelman 0, 146997304618SKris Buschelman 0, 147097304618SKris Buschelman 0, 147197304618SKris Buschelman /*75*/ 0, 147297304618SKris Buschelman 0, 147397304618SKris Buschelman 0, 147497304618SKris Buschelman 0, 147597304618SKris Buschelman 0, 147697304618SKris Buschelman /*80*/ 0, 147797304618SKris Buschelman 0, 147897304618SKris Buschelman 0, 147997304618SKris Buschelman 0, 148097304618SKris Buschelman 0, 148197304618SKris Buschelman /*85*/ MatLoad_SeqDense}; 148290ace30eSBarry Smith 14834a2ae208SSatish Balay #undef __FUNCT__ 14844a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 14854b828684SBarry Smith /*@C 1486fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1487d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1488d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1489289bc588SBarry Smith 1490db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1491db81eaa0SLois Curfman McInnes 149220563c6bSBarry Smith Input Parameters: 1493db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 14940c775827SLois Curfman McInnes . m - number of rows 149518f449edSLois Curfman McInnes . n - number of columns 1496db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1497dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 149820563c6bSBarry Smith 149920563c6bSBarry Smith Output Parameter: 150044cd7ae7SLois Curfman McInnes . A - the matrix 150120563c6bSBarry Smith 1502b259b22eSLois Curfman McInnes Notes: 150318f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 150418f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1505b4fd4287SBarry Smith set data=PETSC_NULL. 150618f449edSLois Curfman McInnes 1507027ccd11SLois Curfman McInnes Level: intermediate 1508027ccd11SLois Curfman McInnes 1509dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1510d65003e9SLois Curfman McInnes 1511db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 151220563c6bSBarry Smith @*/ 151387828ca2SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1514289bc588SBarry Smith { 1515273d9f13SBarry Smith int ierr; 15163b2fbd54SBarry Smith 15173a40ed3dSBarry Smith PetscFunctionBegin; 1518273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1519273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1520273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1521273d9f13SBarry Smith PetscFunctionReturn(0); 1522273d9f13SBarry Smith } 1523273d9f13SBarry Smith 15244a2ae208SSatish Balay #undef __FUNCT__ 15254a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1526273d9f13SBarry Smith /*@C 1527273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1528273d9f13SBarry Smith 1529273d9f13SBarry Smith Collective on MPI_Comm 1530273d9f13SBarry Smith 1531273d9f13SBarry Smith Input Parameters: 1532273d9f13SBarry Smith + A - the matrix 1533273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1534273d9f13SBarry Smith 1535273d9f13SBarry Smith Notes: 1536273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1537273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1538273d9f13SBarry Smith set data=PETSC_NULL. 1539273d9f13SBarry Smith 1540273d9f13SBarry Smith Level: intermediate 1541273d9f13SBarry Smith 1542273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1543273d9f13SBarry Smith 1544273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1545273d9f13SBarry Smith @*/ 1546ca01db9bSBarry Smith int MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1547273d9f13SBarry Smith { 1548ca01db9bSBarry Smith int ierr,(*f)(Mat,PetscScalar[]); 1549a23d5eceSKris Buschelman 1550a23d5eceSKris Buschelman PetscFunctionBegin; 1551a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1552a23d5eceSKris Buschelman if (f) { 1553a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1554a23d5eceSKris Buschelman } 1555a23d5eceSKris Buschelman PetscFunctionReturn(0); 1556a23d5eceSKris Buschelman } 1557a23d5eceSKris Buschelman 1558a23d5eceSKris Buschelman EXTERN_C_BEGIN 1559a23d5eceSKris Buschelman #undef __FUNCT__ 1560a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1561a23d5eceSKris Buschelman int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1562a23d5eceSKris Buschelman { 1563273d9f13SBarry Smith Mat_SeqDense *b; 1564273d9f13SBarry Smith int ierr; 1565273d9f13SBarry Smith 1566273d9f13SBarry Smith PetscFunctionBegin; 1567273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1568273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1569273d9f13SBarry Smith if (!data) { 157087828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 157187828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1572273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 157387828ca2SBarry Smith PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1574273d9f13SBarry Smith } else { /* user-allocated storage */ 1575273d9f13SBarry Smith b->v = data; 1576273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1577273d9f13SBarry Smith } 1578273d9f13SBarry Smith PetscFunctionReturn(0); 1579273d9f13SBarry Smith } 1580a23d5eceSKris Buschelman EXTERN_C_END 1581273d9f13SBarry Smith 15821b807ce4Svictorle #undef __FUNCT__ 15831b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 15841b807ce4Svictorle /*@C 15851b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 15861b807ce4Svictorle 15871b807ce4Svictorle Input parameter: 15881b807ce4Svictorle + A - the matrix 15891b807ce4Svictorle - lda - the leading dimension 15901b807ce4Svictorle 15911b807ce4Svictorle Notes: 15921b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 15931b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 15941b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 15951b807ce4Svictorle 15961b807ce4Svictorle Level: intermediate 15971b807ce4Svictorle 15981b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 15991b807ce4Svictorle 16001b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 16011b807ce4Svictorle @*/ 16021b807ce4Svictorle int MatSeqDenseSetLDA(Mat B,int lda) 16031b807ce4Svictorle { 16041b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 16051b807ce4Svictorle PetscFunctionBegin; 16061b807ce4Svictorle if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension"); 16071b807ce4Svictorle b->lda = lda; 16081b807ce4Svictorle PetscFunctionReturn(0); 16091b807ce4Svictorle } 16101b807ce4Svictorle 16110bad9183SKris Buschelman /*MC 1612fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 16130bad9183SKris Buschelman 16140bad9183SKris Buschelman Options Database Keys: 16150bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 16160bad9183SKris Buschelman 16170bad9183SKris Buschelman Level: beginner 16180bad9183SKris Buschelman 16190bad9183SKris Buschelman .seealso: MatCreateSeqDense 16200bad9183SKris Buschelman M*/ 16210bad9183SKris Buschelman 1622273d9f13SBarry Smith EXTERN_C_BEGIN 16234a2ae208SSatish Balay #undef __FUNCT__ 16244a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1625273d9f13SBarry Smith int MatCreate_SeqDense(Mat B) 1626273d9f13SBarry Smith { 1627273d9f13SBarry Smith Mat_SeqDense *b; 1628273d9f13SBarry Smith int ierr,size; 1629273d9f13SBarry Smith 1630273d9f13SBarry Smith PetscFunctionBegin; 1631273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 163229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 163355659b69SBarry Smith 1634273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1635273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1636273d9f13SBarry Smith 1637b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1638549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1639549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 164044cd7ae7SLois Curfman McInnes B->factor = 0; 164190f02eecSBarry Smith B->mapping = 0; 1642b0a32e0cSBarry Smith PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 164344cd7ae7SLois Curfman McInnes B->data = (void*)b; 164418f449edSLois Curfman McInnes 16458a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 16468a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1647a5ae1ecdSBarry Smith 164844cd7ae7SLois Curfman McInnes b->pivots = 0; 1649273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1650273d9f13SBarry Smith b->v = 0; 16511b807ce4Svictorle b->lda = B->m; 16524e220ebcSLois Curfman McInnes 1653a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1654a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 1655a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 16563a40ed3dSBarry Smith PetscFunctionReturn(0); 1657289bc588SBarry Smith } 1658273d9f13SBarry Smith EXTERN_C_END 1659