167e560aaSBarry Smith /* 267e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 367e560aaSBarry Smith */ 4289bc588SBarry Smith 570f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 6b4c862e0SSatish Balay #include "petscblaslapack.h" 7289bc588SBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 10dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str) 111987afe7SBarry Smith { 121987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 1313f74950SBarry Smith PetscInt j; 144ce68768SBarry Smith PetscBLASInt N = (PetscBLASInt)X->m*X->n,m=(PetscBLASInt)X->m,ldax = x->lda,lday=y->lda,one = 1; 15*efee365bSSatish Balay PetscErrorCode ierr; 163a40ed3dSBarry Smith 173a40ed3dSBarry Smith PetscFunctionBegin; 18a5ce6ee0Svictorle if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 19a5ce6ee0Svictorle if (ldax>m || lday>m) { 20a5ce6ee0Svictorle for (j=0; j<X->n; j++) { 2171044d3cSBarry Smith BLASaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one); 22a5ce6ee0Svictorle } 23a5ce6ee0Svictorle } else { 2471044d3cSBarry Smith BLASaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one); 25a5ce6ee0Svictorle } 26*efee365bSSatish Balay ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr); 273a40ed3dSBarry Smith PetscFunctionReturn(0); 281987afe7SBarry Smith } 291987afe7SBarry Smith 304a2ae208SSatish Balay #undef __FUNCT__ 314a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 32dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 33289bc588SBarry Smith { 344e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3513f74950SBarry Smith PetscInt i,N = A->m*A->n,count = 0; 3687828ca2SBarry Smith PetscScalar *v = a->v; 373a40ed3dSBarry Smith 383a40ed3dSBarry Smith PetscFunctionBegin; 39289bc588SBarry Smith for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 404e220ebcSLois Curfman McInnes 41273d9f13SBarry Smith info->rows_global = (double)A->m; 42273d9f13SBarry Smith info->columns_global = (double)A->n; 43273d9f13SBarry Smith info->rows_local = (double)A->m; 44273d9f13SBarry Smith info->columns_local = (double)A->n; 454e220ebcSLois Curfman McInnes info->block_size = 1.0; 464e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 474e220ebcSLois Curfman McInnes info->nz_used = (double)count; 484e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 494e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 504e220ebcSLois Curfman McInnes info->mallocs = 0; 514e220ebcSLois Curfman McInnes info->memory = A->mem; 524e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 534e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 544e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 554e220ebcSLois Curfman McInnes 563a40ed3dSBarry Smith PetscFunctionReturn(0); 57289bc588SBarry Smith } 58289bc588SBarry Smith 594a2ae208SSatish Balay #undef __FUNCT__ 604a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 61dfbe8321SBarry Smith PetscErrorCode MatScale_SeqDense(const PetscScalar *alpha,Mat A) 6280cd9d93SLois Curfman McInnes { 63273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 644ce68768SBarry Smith PetscBLASInt one = 1,lda = a->lda,j,nz; 65*efee365bSSatish Balay PetscErrorCode ierr; 6680cd9d93SLois Curfman McInnes 673a40ed3dSBarry Smith PetscFunctionBegin; 68a5ce6ee0Svictorle if (lda>A->m) { 694ce68768SBarry Smith nz = (PetscBLASInt)A->m; 70a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 7171044d3cSBarry Smith BLASscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one); 72a5ce6ee0Svictorle } 73a5ce6ee0Svictorle } else { 744ce68768SBarry Smith nz = (PetscBLASInt)A->m*A->n; 7571044d3cSBarry Smith BLASscal_(&nz,(PetscScalar*)alpha,a->v,&one); 76a5ce6ee0Svictorle } 77*efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 783a40ed3dSBarry Smith PetscFunctionReturn(0); 7980cd9d93SLois Curfman McInnes } 8080cd9d93SLois Curfman McInnes 81289bc588SBarry Smith /* ---------------------------------------------------------------*/ 826831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 836831982aSBarry Smith rather than put it in the Mat->row slot.*/ 844a2ae208SSatish Balay #undef __FUNCT__ 854a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 86dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 87289bc588SBarry Smith { 88a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 8981f479a6Svictorle PetscFunctionBegin; 90a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 91a07ab388SBarry Smith #else 92c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 936849ba73SBarry Smith PetscErrorCode ierr; 944ce68768SBarry Smith PetscBLASInt n = (PetscBLASInt)A->n,m = (PetscBLASInt)A->m,info; 9555659b69SBarry Smith 963a40ed3dSBarry Smith PetscFunctionBegin; 97289bc588SBarry Smith if (!mat->pivots) { 984ce68768SBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 9952e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,A->m*sizeof(PetscBLASInt));CHKERRQ(ierr); 100289bc588SBarry Smith } 101f1af5d2fSBarry Smith A->factor = FACTOR_LU; 102273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 10371044d3cSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 10429bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 10529bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 106*efee365bSSatish Balay ierr = PetscLogFlops((2*A->n*A->n*A->n)/3);CHKERRQ(ierr); 107a07ab388SBarry Smith #endif 1083a40ed3dSBarry Smith PetscFunctionReturn(0); 109289bc588SBarry Smith } 1106ee01492SSatish Balay 1114a2ae208SSatish Balay #undef __FUNCT__ 1124a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 113dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 11402cad45dSBarry Smith { 11502cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 1166849ba73SBarry Smith PetscErrorCode ierr; 11713f74950SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 11802cad45dSBarry Smith Mat newi; 11902cad45dSBarry Smith 1203a40ed3dSBarry Smith PetscFunctionBegin; 1215c5985e7SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&newi);CHKERRQ(ierr); 1225c5985e7SKris Buschelman ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr); 1235c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 1245609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 125a5ce6ee0Svictorle l = (Mat_SeqDense*)newi->data; 126a5ce6ee0Svictorle if (lda>A->m) { 127a5ce6ee0Svictorle m = A->m; 128a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 129a5ce6ee0Svictorle ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 130a5ce6ee0Svictorle } 131a5ce6ee0Svictorle } else { 13287828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 13302cad45dSBarry Smith } 134a5ce6ee0Svictorle } 13502cad45dSBarry Smith newi->assembled = PETSC_TRUE; 13602cad45dSBarry Smith *newmat = newi; 1373a40ed3dSBarry Smith PetscFunctionReturn(0); 13802cad45dSBarry Smith } 13902cad45dSBarry Smith 1404a2ae208SSatish Balay #undef __FUNCT__ 1414a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 142dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 143289bc588SBarry Smith { 144dfbe8321SBarry Smith PetscErrorCode ierr; 1453a40ed3dSBarry Smith 1463a40ed3dSBarry Smith PetscFunctionBegin; 1475609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1483a40ed3dSBarry Smith PetscFunctionReturn(0); 149289bc588SBarry Smith } 1506ee01492SSatish Balay 1514a2ae208SSatish Balay #undef __FUNCT__ 1524a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 153af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 154289bc588SBarry Smith { 15502cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 1566849ba73SBarry Smith PetscErrorCode ierr; 15713f74950SBarry Smith PetscInt lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j; 1584482741eSBarry Smith MatFactorInfo info; 1593a40ed3dSBarry Smith 1603a40ed3dSBarry Smith PetscFunctionBegin; 16102cad45dSBarry Smith /* copy the numerical values */ 1621b807ce4Svictorle if (lda1>m || lda2>m ) { 1631b807ce4Svictorle for (j=0; j<n; j++) { 1641b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1651b807ce4Svictorle } 1661b807ce4Svictorle } else { 16787828ca2SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1681b807ce4Svictorle } 16902cad45dSBarry Smith (*fact)->factor = 0; 1704482741eSBarry Smith ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr); 1713a40ed3dSBarry Smith PetscFunctionReturn(0); 172289bc588SBarry Smith } 1736ee01492SSatish Balay 1744a2ae208SSatish Balay #undef __FUNCT__ 1754a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 176dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 177289bc588SBarry Smith { 178dfbe8321SBarry Smith PetscErrorCode ierr; 1793a40ed3dSBarry Smith 1803a40ed3dSBarry Smith PetscFunctionBegin; 181ceb03754SKris Buschelman ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr); 1823a40ed3dSBarry Smith PetscFunctionReturn(0); 183289bc588SBarry Smith } 1846ee01492SSatish Balay 1854a2ae208SSatish Balay #undef __FUNCT__ 1864a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 187dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 188289bc588SBarry Smith { 189a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 190a07ab388SBarry Smith PetscFunctionBegin; 191a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 192a07ab388SBarry Smith #else 193c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1946849ba73SBarry Smith PetscErrorCode ierr; 1954ce68768SBarry Smith PetscBLASInt n = (PetscBLASInt)A->n,info; 19655659b69SBarry Smith 1973a40ed3dSBarry Smith PetscFunctionBegin; 198752f5567SLois Curfman McInnes if (mat->pivots) { 199606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 20052e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-A->m*sizeof(PetscInt));CHKERRQ(ierr); 201752f5567SLois Curfman McInnes mat->pivots = 0; 202752f5567SLois Curfman McInnes } 203273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 20471044d3cSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 20577431f27SBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 206c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 207*efee365bSSatish Balay ierr = PetscLogFlops((A->n*A->n*A->n)/3);CHKERRQ(ierr); 208a07ab388SBarry Smith #endif 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 210289bc588SBarry Smith } 211289bc588SBarry Smith 2124a2ae208SSatish Balay #undef __FUNCT__ 2130b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 214af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 2150b4b3355SBarry Smith { 216dfbe8321SBarry Smith PetscErrorCode ierr; 21715e8a5b3SHong Zhang MatFactorInfo info; 2180b4b3355SBarry Smith 2190b4b3355SBarry Smith PetscFunctionBegin; 22015e8a5b3SHong Zhang info.fill = 1.0; 22115e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 2220b4b3355SBarry Smith PetscFunctionReturn(0); 2230b4b3355SBarry Smith } 2240b4b3355SBarry Smith 2250b4b3355SBarry Smith #undef __FUNCT__ 2264a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 227dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 228289bc588SBarry Smith { 229c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2306849ba73SBarry Smith PetscErrorCode ierr; 2314ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m, one = 1,info; 23287828ca2SBarry Smith PetscScalar *x,*y; 23367e560aaSBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 235273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2361ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2371ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 23887828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 239c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 240ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 24180444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 242ae7cfcebSSatish Balay #else 24371044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2444ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve"); 245ae7cfcebSSatish Balay #endif 2467a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 247ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 24880444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 249ae7cfcebSSatish Balay #else 25071044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 2514ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve"); 252ae7cfcebSSatish Balay #endif 253289bc588SBarry Smith } 25429bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2551ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2561ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 257*efee365bSSatish Balay ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr); 2583a40ed3dSBarry Smith PetscFunctionReturn(0); 259289bc588SBarry Smith } 2606ee01492SSatish Balay 2614a2ae208SSatish Balay #undef __FUNCT__ 2624a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 263dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 264da3a660dSBarry Smith { 265c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 266dfbe8321SBarry Smith PetscErrorCode ierr; 2674ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt) A->m,one = 1,info; 26887828ca2SBarry Smith PetscScalar *x,*y; 26967e560aaSBarry Smith 2703a40ed3dSBarry Smith PetscFunctionBegin; 271273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 2721ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2731ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 27487828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 275752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 276da3a660dSBarry Smith if (mat->pivots) { 277ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 27880444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 279ae7cfcebSSatish Balay #else 28071044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2814ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 282ae7cfcebSSatish Balay #endif 2837a97a34bSBarry Smith } else { 284ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 28580444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 286ae7cfcebSSatish Balay #else 28771044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 2884ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 289ae7cfcebSSatish Balay #endif 290da3a660dSBarry Smith } 2911ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2921ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 293*efee365bSSatish Balay ierr = PetscLogFlops(2*A->n*A->n - A->n);CHKERRQ(ierr); 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 295da3a660dSBarry Smith } 2966ee01492SSatish Balay 2974a2ae208SSatish Balay #undef __FUNCT__ 2984a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 299dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 300da3a660dSBarry Smith { 301c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 302dfbe8321SBarry Smith PetscErrorCode ierr; 3034ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m,one = 1,info; 30487828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 305da3a660dSBarry Smith Vec tmp = 0; 30667e560aaSBarry Smith 3073a40ed3dSBarry Smith PetscFunctionBegin; 3081ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3091ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 310273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 311da3a660dSBarry Smith if (yy == zz) { 31278b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 31352e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 31478b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 315da3a660dSBarry Smith } 31687828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 317752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 318da3a660dSBarry Smith if (mat->pivots) { 319ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 32080444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 321ae7cfcebSSatish Balay #else 32271044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3232ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 324ae7cfcebSSatish Balay #endif 325a8c6a408SBarry Smith } else { 326ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 32780444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 328ae7cfcebSSatish Balay #else 32971044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3302ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 331ae7cfcebSSatish Balay #endif 332da3a660dSBarry Smith } 333a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 334a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 3351ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3361ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 337*efee365bSSatish Balay ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr); 3383a40ed3dSBarry Smith PetscFunctionReturn(0); 339da3a660dSBarry Smith } 34067e560aaSBarry Smith 3414a2ae208SSatish Balay #undef __FUNCT__ 3424a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 343dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 344da3a660dSBarry Smith { 345c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3466849ba73SBarry Smith PetscErrorCode ierr; 3474ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m,one = 1,info; 34887828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 349da3a660dSBarry Smith Vec tmp; 35067e560aaSBarry Smith 3513a40ed3dSBarry Smith PetscFunctionBegin; 352273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 3531ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3541ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 355da3a660dSBarry Smith if (yy == zz) { 35678b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 35752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 35878b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 359da3a660dSBarry Smith } 36087828ca2SBarry Smith ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 361752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 362da3a660dSBarry Smith if (mat->pivots) { 363ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 36480444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 365ae7cfcebSSatish Balay #else 36671044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3672ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 368ae7cfcebSSatish Balay #endif 3693a40ed3dSBarry Smith } else { 370ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 37180444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 372ae7cfcebSSatish Balay #else 37371044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3742ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 375ae7cfcebSSatish Balay #endif 376da3a660dSBarry Smith } 37790f02eecSBarry Smith if (tmp) { 37890f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 37990f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3803a40ed3dSBarry Smith } else { 38190f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 38290f02eecSBarry Smith } 3831ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3841ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 385*efee365bSSatish Balay ierr = PetscLogFlops(2*A->n*A->n);CHKERRQ(ierr); 3863a40ed3dSBarry Smith PetscFunctionReturn(0); 387da3a660dSBarry Smith } 388289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3894a2ae208SSatish Balay #undef __FUNCT__ 3904a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 39113f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 392289bc588SBarry Smith { 393c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 39487828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 395dfbe8321SBarry Smith PetscErrorCode ierr; 39613f74950SBarry Smith PetscInt m = A->m,i; 397aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3984ce68768SBarry Smith PetscBLASInt bm = (PetscBLASInt)m, o = 1; 399bc1b551cSSatish Balay #endif 400289bc588SBarry Smith 4013a40ed3dSBarry Smith PetscFunctionBegin; 402289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 40371044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 404bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 405289bc588SBarry Smith } 4061ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4071ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 408b965ef7fSBarry Smith its = its*lits; 40977431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 410289bc588SBarry Smith while (its--) { 411289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 412289bc588SBarry Smith for (i=0; i<m; i++) { 413aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 414f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 415f1747703SBarry Smith not happy about returning a double complex */ 41613f74950SBarry Smith PetscInt _i; 41787828ca2SBarry Smith PetscScalar sum = b[i]; 418f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4193f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 420f1747703SBarry Smith } 421f1747703SBarry Smith xt = sum; 422f1747703SBarry Smith #else 42371044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 424f1747703SBarry Smith #endif 42555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 426289bc588SBarry Smith } 427289bc588SBarry Smith } 428289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 429289bc588SBarry Smith for (i=m-1; i>=0; i--) { 430aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 431f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 432f1747703SBarry Smith not happy about returning a double complex */ 43313f74950SBarry Smith PetscInt _i; 43487828ca2SBarry Smith PetscScalar sum = b[i]; 435f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4363f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 437f1747703SBarry Smith } 438f1747703SBarry Smith xt = sum; 439f1747703SBarry Smith #else 44071044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 441f1747703SBarry Smith #endif 44255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 443289bc588SBarry Smith } 444289bc588SBarry Smith } 445289bc588SBarry Smith } 4461ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4471ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4483a40ed3dSBarry Smith PetscFunctionReturn(0); 449289bc588SBarry Smith } 450289bc588SBarry Smith 451289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4524a2ae208SSatish Balay #undef __FUNCT__ 4534a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 454dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 455289bc588SBarry Smith { 456c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 45787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 458dfbe8321SBarry Smith PetscErrorCode ierr; 4594ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n,_One=1; 460ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4613a40ed3dSBarry Smith 4623a40ed3dSBarry Smith PetscFunctionBegin; 463273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4641ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4651ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 46671044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 4671ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4681ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 469*efee365bSSatish Balay ierr = PetscLogFlops(2*A->m*A->n - A->n);CHKERRQ(ierr); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 471289bc588SBarry Smith } 4726ee01492SSatish Balay 4734a2ae208SSatish Balay #undef __FUNCT__ 4744a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 475dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 476289bc588SBarry Smith { 477c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 47887828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 479dfbe8321SBarry Smith PetscErrorCode ierr; 4804ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1; 4813a40ed3dSBarry Smith 4823a40ed3dSBarry Smith PetscFunctionBegin; 483273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 4841ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4851ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 48671044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 4871ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4881ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 489*efee365bSSatish Balay ierr = PetscLogFlops(2*A->m*A->n - A->m);CHKERRQ(ierr); 4903a40ed3dSBarry Smith PetscFunctionReturn(0); 491289bc588SBarry Smith } 4926ee01492SSatish Balay 4934a2ae208SSatish Balay #undef __FUNCT__ 4944a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 495dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 496289bc588SBarry Smith { 497c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 49887828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 499dfbe8321SBarry Smith PetscErrorCode ierr; 5004ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1; 5013a40ed3dSBarry Smith 5023a40ed3dSBarry Smith PetscFunctionBegin; 503273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 504600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5051ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5061ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 50771044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5081ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5091ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 510*efee365bSSatish Balay ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr); 5113a40ed3dSBarry Smith PetscFunctionReturn(0); 512289bc588SBarry Smith } 5136ee01492SSatish Balay 5144a2ae208SSatish Balay #undef __FUNCT__ 5154a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 516dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 517289bc588SBarry Smith { 518c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 51987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 520dfbe8321SBarry Smith PetscErrorCode ierr; 5214ce68768SBarry Smith PetscBLASInt m = (PetscBLASInt)A->m, n = (PetscBLASInt)A->n, _One=1; 52287828ca2SBarry Smith PetscScalar _DOne=1.0; 5233a40ed3dSBarry Smith 5243a40ed3dSBarry Smith PetscFunctionBegin; 525273d9f13SBarry Smith if (!A->m || !A->n) PetscFunctionReturn(0); 526600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5271ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5281ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 52971044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5301ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5311ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 532*efee365bSSatish Balay ierr = PetscLogFlops(2*A->m*A->n);CHKERRQ(ierr); 5333a40ed3dSBarry Smith PetscFunctionReturn(0); 534289bc588SBarry Smith } 535289bc588SBarry Smith 536289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5374a2ae208SSatish Balay #undef __FUNCT__ 5384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 53913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 540289bc588SBarry Smith { 541c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54287828ca2SBarry Smith PetscScalar *v; 5436849ba73SBarry Smith PetscErrorCode ierr; 54413f74950SBarry Smith PetscInt i; 54567e560aaSBarry Smith 5463a40ed3dSBarry Smith PetscFunctionBegin; 547273d9f13SBarry Smith *ncols = A->n; 548289bc588SBarry Smith if (cols) { 54913f74950SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 550273d9f13SBarry Smith for (i=0; i<A->n; i++) (*cols)[i] = i; 551289bc588SBarry Smith } 552289bc588SBarry Smith if (vals) { 55387828ca2SBarry Smith ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 554289bc588SBarry Smith v = mat->v + row; 5551b807ce4Svictorle for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;} 556289bc588SBarry Smith } 5573a40ed3dSBarry Smith PetscFunctionReturn(0); 558289bc588SBarry Smith } 5596ee01492SSatish Balay 5604a2ae208SSatish Balay #undef __FUNCT__ 5614a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 56213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 563289bc588SBarry Smith { 564dfbe8321SBarry Smith PetscErrorCode ierr; 565606d414cSSatish Balay PetscFunctionBegin; 566606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 567606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 569289bc588SBarry Smith } 570289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5714a2ae208SSatish Balay #undef __FUNCT__ 5724a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 57313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 574289bc588SBarry Smith { 575c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57613f74950SBarry Smith PetscInt i,j,idx=0; 577d6dfbf8fSBarry Smith 5783a40ed3dSBarry Smith PetscFunctionBegin; 579289bc588SBarry Smith if (!mat->roworiented) { 580dbb450caSBarry Smith if (addv == INSERT_VALUES) { 581289bc588SBarry Smith for (j=0; j<n; j++) { 582cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 5832515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 58477431f27SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1); 58558804f6dSBarry Smith #endif 586289bc588SBarry Smith for (i=0; i<m; i++) { 587cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 5882515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 58977431f27SBarry Smith if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1); 59058804f6dSBarry Smith #endif 591cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 592289bc588SBarry Smith } 593289bc588SBarry Smith } 5943a40ed3dSBarry Smith } else { 595289bc588SBarry Smith for (j=0; j<n; j++) { 596cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 5972515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 59877431f27SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1); 59958804f6dSBarry Smith #endif 600289bc588SBarry Smith for (i=0; i<m; i++) { 601cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6022515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 60377431f27SBarry Smith if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1); 60458804f6dSBarry Smith #endif 605cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 606289bc588SBarry Smith } 607289bc588SBarry Smith } 608289bc588SBarry Smith } 6093a40ed3dSBarry Smith } else { 610dbb450caSBarry Smith if (addv == INSERT_VALUES) { 611e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 612cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6132515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61477431f27SBarry Smith if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1); 61558804f6dSBarry Smith #endif 616e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 617cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6182515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61977431f27SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1); 62058804f6dSBarry Smith #endif 621cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 622e8d4e0b9SBarry Smith } 623e8d4e0b9SBarry Smith } 6243a40ed3dSBarry Smith } else { 625289bc588SBarry Smith for (i=0; i<m; i++) { 626cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6272515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 62877431f27SBarry Smith if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->m-1); 62958804f6dSBarry Smith #endif 630289bc588SBarry Smith for (j=0; j<n; j++) { 631cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6322515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 63377431f27SBarry Smith if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->n-1); 63458804f6dSBarry Smith #endif 635cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 636289bc588SBarry Smith } 637289bc588SBarry Smith } 638289bc588SBarry Smith } 639e8d4e0b9SBarry Smith } 6403a40ed3dSBarry Smith PetscFunctionReturn(0); 641289bc588SBarry Smith } 642e8d4e0b9SBarry Smith 6434a2ae208SSatish Balay #undef __FUNCT__ 6444a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 64513f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 646ae80bb75SLois Curfman McInnes { 647ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64813f74950SBarry Smith PetscInt i,j; 64987828ca2SBarry Smith PetscScalar *vpt = v; 650ae80bb75SLois Curfman McInnes 6513a40ed3dSBarry Smith PetscFunctionBegin; 652ae80bb75SLois Curfman McInnes /* row-oriented output */ 653ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 654ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 6551b807ce4Svictorle *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 656ae80bb75SLois Curfman McInnes } 657ae80bb75SLois Curfman McInnes } 6583a40ed3dSBarry Smith PetscFunctionReturn(0); 659ae80bb75SLois Curfman McInnes } 660ae80bb75SLois Curfman McInnes 661289bc588SBarry Smith /* -----------------------------------------------------------------*/ 662289bc588SBarry Smith 663e090d566SSatish Balay #include "petscsys.h" 66488e32edaSLois Curfman McInnes 6654a2ae208SSatish Balay #undef __FUNCT__ 6664a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 667dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A) 66888e32edaSLois Curfman McInnes { 66988e32edaSLois Curfman McInnes Mat_SeqDense *a; 67088e32edaSLois Curfman McInnes Mat B; 6716849ba73SBarry Smith PetscErrorCode ierr; 67213f74950SBarry Smith PetscInt *scols,i,j,nz,header[4]; 67313f74950SBarry Smith int fd; 67413f74950SBarry Smith PetscMPIInt size; 67513f74950SBarry Smith PetscInt *rowlengths = 0,M,N,*cols; 67687828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 67719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 67888e32edaSLois Curfman McInnes 6793a40ed3dSBarry Smith PetscFunctionBegin; 680d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 68129bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 682b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 6830752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 684552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 68588e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 68688e32edaSLois Curfman McInnes 68790ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 6885c5985e7SKris Buschelman ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr); 689be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 6905c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 69190ace30eSBarry Smith B = *A; 69290ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 6934a41a4d3SSatish Balay v = a->v; 6944a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 6954a41a4d3SSatish Balay from row major to column major */ 696b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 69790ace30eSBarry Smith /* read in nonzero values */ 6984a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6994a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 7004a41a4d3SSatish Balay for (j=0; j<N; j++) { 7014a41a4d3SSatish Balay for (i=0; i<M; i++) { 7024a41a4d3SSatish Balay *v++ =w[i*N+j]; 7034a41a4d3SSatish Balay } 7044a41a4d3SSatish Balay } 705606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 7066d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7076d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70890ace30eSBarry Smith } else { 70988e32edaSLois Curfman McInnes /* read row lengths */ 71013f74950SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 7110752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 71288e32edaSLois Curfman McInnes 71388e32edaSLois Curfman McInnes /* create our matrix */ 7145c5985e7SKris Buschelman ierr = MatCreate(comm,M,N,M,N,A);CHKERRQ(ierr); 715be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7165c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 71788e32edaSLois Curfman McInnes B = *A; 71888e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 71988e32edaSLois Curfman McInnes v = a->v; 72088e32edaSLois Curfman McInnes 72188e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 72213f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 723b0a32e0cSBarry Smith cols = scols; 7240752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 72587828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 726b0a32e0cSBarry Smith vals = svals; 7270752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 72888e32edaSLois Curfman McInnes 72988e32edaSLois Curfman McInnes /* insert into matrix */ 73088e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 73188e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 73288e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 73388e32edaSLois Curfman McInnes } 734606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 735606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 736606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 73788e32edaSLois Curfman McInnes 7386d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7396d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74090ace30eSBarry Smith } 7413a40ed3dSBarry Smith PetscFunctionReturn(0); 74288e32edaSLois Curfman McInnes } 74388e32edaSLois Curfman McInnes 744e090d566SSatish Balay #include "petscsys.h" 745932b0c3eSLois Curfman McInnes 7464a2ae208SSatish Balay #undef __FUNCT__ 7474a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 7486849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 749289bc588SBarry Smith { 750932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 751dfbe8321SBarry Smith PetscErrorCode ierr; 75213f74950SBarry Smith PetscInt i,j; 753fb9695e5SSatish Balay char *name; 75487828ca2SBarry Smith PetscScalar *v; 755f3ef73ceSBarry Smith PetscViewerFormat format; 756932b0c3eSLois Curfman McInnes 7573a40ed3dSBarry Smith PetscFunctionBegin; 758435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 759b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 760456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7613a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 762fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 763b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 764273d9f13SBarry Smith for (i=0; i<A->m; i++) { 76544cd7ae7SLois Curfman McInnes v = a->v + i; 76677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 767273d9f13SBarry Smith for (j=0; j<A->n; j++) { 768aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 769329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 77077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 771329f5518SBarry Smith } else if (PetscRealPart(*v)) { 77277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7736831982aSBarry Smith } 77480cd9d93SLois Curfman McInnes #else 7756831982aSBarry Smith if (*v) { 77677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,*v);CHKERRQ(ierr); 7776831982aSBarry Smith } 77880cd9d93SLois Curfman McInnes #endif 7791b807ce4Svictorle v += a->lda; 78080cd9d93SLois Curfman McInnes } 781b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 78280cd9d93SLois Curfman McInnes } 783b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 7843a40ed3dSBarry Smith } else { 785b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 786aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 787ffac6cdbSBarry Smith PetscTruth allreal = PETSC_TRUE; 78847989497SBarry Smith /* determine if matrix has all real values */ 78947989497SBarry Smith v = a->v; 790201f5f94SSatish Balay for (i=0; i<A->m*A->n; i++) { 791ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 79247989497SBarry Smith } 79347989497SBarry Smith #endif 794fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 7953a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 79677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->m,A->n);CHKERRQ(ierr); 79777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->m,A->n);CHKERRQ(ierr); 798fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 799ffac6cdbSBarry Smith } 800ffac6cdbSBarry Smith 801273d9f13SBarry Smith for (i=0; i<A->m; i++) { 802932b0c3eSLois Curfman McInnes v = a->v + i; 803273d9f13SBarry Smith for (j=0; j<A->n; j++) { 804aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 80547989497SBarry Smith if (allreal) { 8061b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 80747989497SBarry Smith } else { 8081b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 80947989497SBarry Smith } 810289bc588SBarry Smith #else 8111b807ce4Svictorle ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 812289bc588SBarry Smith #endif 8131b807ce4Svictorle v += a->lda; 814289bc588SBarry Smith } 815b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 816289bc588SBarry Smith } 817fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 818b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 819ffac6cdbSBarry Smith } 820b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 821da3a660dSBarry Smith } 822b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8233a40ed3dSBarry Smith PetscFunctionReturn(0); 824289bc588SBarry Smith } 825289bc588SBarry Smith 8264a2ae208SSatish Balay #undef __FUNCT__ 8274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 8286849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 829932b0c3eSLois Curfman McInnes { 830932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 8316849ba73SBarry Smith PetscErrorCode ierr; 83213f74950SBarry Smith int fd; 83313f74950SBarry Smith PetscInt ict,j,n = A->n,m = A->m,i,*col_lens,nz = m*n; 83487828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 835f3ef73ceSBarry Smith PetscViewerFormat format; 836932b0c3eSLois Curfman McInnes 8373a40ed3dSBarry Smith PetscFunctionBegin; 838b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 83990ace30eSBarry Smith 840b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 841fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 84290ace30eSBarry Smith /* store the matrix as a dense matrix */ 84313f74950SBarry Smith ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 844552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 84590ace30eSBarry Smith col_lens[1] = m; 84690ace30eSBarry Smith col_lens[2] = n; 84790ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8486f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 849606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 85090ace30eSBarry Smith 85190ace30eSBarry Smith /* write out matrix, by rows */ 85287828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 85390ace30eSBarry Smith v = a->v; 85490ace30eSBarry Smith for (i=0; i<m; i++) { 85590ace30eSBarry Smith for (j=0; j<n; j++) { 85690ace30eSBarry Smith vals[i + j*m] = *v++; 85790ace30eSBarry Smith } 85890ace30eSBarry Smith } 8596f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 860606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 86190ace30eSBarry Smith } else { 86213f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 863552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 864932b0c3eSLois Curfman McInnes col_lens[1] = m; 865932b0c3eSLois Curfman McInnes col_lens[2] = n; 866932b0c3eSLois Curfman McInnes col_lens[3] = nz; 867932b0c3eSLois Curfman McInnes 868932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 869932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8706f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 871932b0c3eSLois Curfman McInnes 872932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 873932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 874932b0c3eSLois Curfman McInnes ict = 0; 875932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 876932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 877932b0c3eSLois Curfman McInnes } 8786f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 879606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 880932b0c3eSLois Curfman McInnes 881932b0c3eSLois Curfman McInnes /* store nonzero values */ 88287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 883932b0c3eSLois Curfman McInnes ict = 0; 884932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 885932b0c3eSLois Curfman McInnes v = a->v + i; 886932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 8871b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 888932b0c3eSLois Curfman McInnes } 889932b0c3eSLois Curfman McInnes } 8906f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 891606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 89290ace30eSBarry Smith } 8933a40ed3dSBarry Smith PetscFunctionReturn(0); 894932b0c3eSLois Curfman McInnes } 895932b0c3eSLois Curfman McInnes 8964a2ae208SSatish Balay #undef __FUNCT__ 8974a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 898dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 899f1af5d2fSBarry Smith { 900f1af5d2fSBarry Smith Mat A = (Mat) Aa; 901f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9026849ba73SBarry Smith PetscErrorCode ierr; 90313f74950SBarry Smith PetscInt m = A->m,n = A->n,color,i,j; 90487828ca2SBarry Smith PetscScalar *v = a->v; 905b0a32e0cSBarry Smith PetscViewer viewer; 906b0a32e0cSBarry Smith PetscDraw popup; 907329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 908f3ef73ceSBarry Smith PetscViewerFormat format; 909f1af5d2fSBarry Smith 910f1af5d2fSBarry Smith PetscFunctionBegin; 911f1af5d2fSBarry Smith 912f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 913b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 914b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 915f1af5d2fSBarry Smith 916f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 917fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 918f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 919b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 920f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 921f1af5d2fSBarry Smith x_l = j; 922f1af5d2fSBarry Smith x_r = x_l + 1.0; 923f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 924f1af5d2fSBarry Smith y_l = m - i - 1.0; 925f1af5d2fSBarry Smith y_r = y_l + 1.0; 926f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 927329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 928b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 929329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 930b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 931f1af5d2fSBarry Smith } else { 932f1af5d2fSBarry Smith continue; 933f1af5d2fSBarry Smith } 934f1af5d2fSBarry Smith #else 935f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 936b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 937f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 938b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 939f1af5d2fSBarry Smith } else { 940f1af5d2fSBarry Smith continue; 941f1af5d2fSBarry Smith } 942f1af5d2fSBarry Smith #endif 943b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 944f1af5d2fSBarry Smith } 945f1af5d2fSBarry Smith } 946f1af5d2fSBarry Smith } else { 947f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 948f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 949f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 950f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 951f1af5d2fSBarry Smith } 952b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 953b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 954b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 955f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 956f1af5d2fSBarry Smith x_l = j; 957f1af5d2fSBarry Smith x_r = x_l + 1.0; 958f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 959f1af5d2fSBarry Smith y_l = m - i - 1.0; 960f1af5d2fSBarry Smith y_r = y_l + 1.0; 961b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 962b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 963f1af5d2fSBarry Smith } 964f1af5d2fSBarry Smith } 965f1af5d2fSBarry Smith } 966f1af5d2fSBarry Smith PetscFunctionReturn(0); 967f1af5d2fSBarry Smith } 968f1af5d2fSBarry Smith 9694a2ae208SSatish Balay #undef __FUNCT__ 9704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 971dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 972f1af5d2fSBarry Smith { 973b0a32e0cSBarry Smith PetscDraw draw; 974f1af5d2fSBarry Smith PetscTruth isnull; 975329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 976dfbe8321SBarry Smith PetscErrorCode ierr; 977f1af5d2fSBarry Smith 978f1af5d2fSBarry Smith PetscFunctionBegin; 979b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 980b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 981abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 982f1af5d2fSBarry Smith 983f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 984273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 985f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 986b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 987b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 988f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 989f1af5d2fSBarry Smith PetscFunctionReturn(0); 990f1af5d2fSBarry Smith } 991f1af5d2fSBarry Smith 9924a2ae208SSatish Balay #undef __FUNCT__ 9934a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 994dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 995932b0c3eSLois Curfman McInnes { 996932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 997dfbe8321SBarry Smith PetscErrorCode ierr; 99832077d6dSBarry Smith PetscTruth issocket,iascii,isbinary,isdraw; 999932b0c3eSLois Curfman McInnes 10003a40ed3dSBarry Smith PetscFunctionBegin; 1001b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 100232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1003fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1004fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 10050f5bd95cSBarry Smith 10060f5bd95cSBarry Smith if (issocket) { 1007634064b4SBarry Smith if (a->lda>A->m) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA"); 1008b0a32e0cSBarry Smith ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 100932077d6dSBarry Smith } else if (iascii) { 10103a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10110f5bd95cSBarry Smith } else if (isbinary) { 10123a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1013f1af5d2fSBarry Smith } else if (isdraw) { 1014f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10155cd90555SBarry Smith } else { 1016958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1017932b0c3eSLois Curfman McInnes } 10183a40ed3dSBarry Smith PetscFunctionReturn(0); 1019932b0c3eSLois Curfman McInnes } 1020289bc588SBarry Smith 10214a2ae208SSatish Balay #undef __FUNCT__ 10224a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1023dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1024289bc588SBarry Smith { 1025ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1026dfbe8321SBarry Smith PetscErrorCode ierr; 102790f02eecSBarry Smith 10283a40ed3dSBarry Smith PetscFunctionBegin; 1029aa482453SBarry Smith #if defined(PETSC_USE_LOG) 103077431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->m,mat->n); 1031a5a9c739SBarry Smith #endif 1032606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 1033606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1034606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1035901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 10363a40ed3dSBarry Smith PetscFunctionReturn(0); 1037289bc588SBarry Smith } 1038289bc588SBarry Smith 10394a2ae208SSatish Balay #undef __FUNCT__ 10404a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1041dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout) 1042289bc588SBarry Smith { 1043c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10446849ba73SBarry Smith PetscErrorCode ierr; 104513f74950SBarry Smith PetscInt k,j,m,n,M; 104687828ca2SBarry Smith PetscScalar *v,tmp; 104748b35521SBarry Smith 10483a40ed3dSBarry Smith PetscFunctionBegin; 10491b807ce4Svictorle v = mat->v; m = A->m; M = mat->lda; n = A->n; 10507c922b88SBarry Smith if (!matout) { /* in place transpose */ 1051a5ce6ee0Svictorle if (m != n) { 1052634064b4SBarry Smith SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1053d64ed03dSBarry Smith } else { 1054d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1055289bc588SBarry Smith for (k=0; k<j; k++) { 10561b807ce4Svictorle tmp = v[j + k*M]; 10571b807ce4Svictorle v[j + k*M] = v[k + j*M]; 10581b807ce4Svictorle v[k + j*M] = tmp; 1059289bc588SBarry Smith } 1060289bc588SBarry Smith } 1061d64ed03dSBarry Smith } 10623a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1063d3e5ee88SLois Curfman McInnes Mat tmat; 1064ec8511deSBarry Smith Mat_SeqDense *tmatd; 106587828ca2SBarry Smith PetscScalar *v2; 1066ea709b57SSatish Balay 10675c5985e7SKris Buschelman ierr = MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);CHKERRQ(ierr); 10685c5985e7SKris Buschelman ierr = MatSetType(tmat,A->type_name);CHKERRQ(ierr); 10695c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1070ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10710de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1072d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10731b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1074d3e5ee88SLois Curfman McInnes } 10756d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10766d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1077d3e5ee88SLois Curfman McInnes *matout = tmat; 107848b35521SBarry Smith } 10793a40ed3dSBarry Smith PetscFunctionReturn(0); 1080289bc588SBarry Smith } 1081289bc588SBarry Smith 10824a2ae208SSatish Balay #undef __FUNCT__ 10834a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1084dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1085289bc588SBarry Smith { 1086c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1087c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 108813f74950SBarry Smith PetscInt i,j; 108987828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 10909ea5d5aeSSatish Balay 10913a40ed3dSBarry Smith PetscFunctionBegin; 1092273d9f13SBarry Smith if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1093273d9f13SBarry Smith if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10941b807ce4Svictorle for (i=0; i<A1->m; i++) { 10951b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 10961b807ce4Svictorle for (j=0; j<A1->n; j++) { 10973a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 10981b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 10991b807ce4Svictorle } 1100289bc588SBarry Smith } 110177c4ece6SBarry Smith *flg = PETSC_TRUE; 11023a40ed3dSBarry Smith PetscFunctionReturn(0); 1103289bc588SBarry Smith } 1104289bc588SBarry Smith 11054a2ae208SSatish Balay #undef __FUNCT__ 11064a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1107dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1108289bc588SBarry Smith { 1109c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1110dfbe8321SBarry Smith PetscErrorCode ierr; 111113f74950SBarry Smith PetscInt i,n,len; 111287828ca2SBarry Smith PetscScalar *x,zero = 0.0; 111344cd7ae7SLois Curfman McInnes 11143a40ed3dSBarry Smith PetscFunctionBegin; 11157a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 11167a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 11171ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1118273d9f13SBarry Smith len = PetscMin(A->m,A->n); 1119273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 112044cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 11211b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1122289bc588SBarry Smith } 11231ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 11243a40ed3dSBarry Smith PetscFunctionReturn(0); 1125289bc588SBarry Smith } 1126289bc588SBarry Smith 11274a2ae208SSatish Balay #undef __FUNCT__ 11284a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1129dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1130289bc588SBarry Smith { 1131c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 113287828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1133dfbe8321SBarry Smith PetscErrorCode ierr; 113413f74950SBarry Smith PetscInt i,j,m = A->m,n = A->n; 113555659b69SBarry Smith 11363a40ed3dSBarry Smith PetscFunctionBegin; 113728988994SBarry Smith if (ll) { 11387a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 11391ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1140273d9f13SBarry Smith if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1141da3a660dSBarry Smith for (i=0; i<m; i++) { 1142da3a660dSBarry Smith x = l[i]; 1143da3a660dSBarry Smith v = mat->v + i; 1144da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1145da3a660dSBarry Smith } 11461ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1147*efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1148da3a660dSBarry Smith } 114928988994SBarry Smith if (rr) { 11507a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 11511ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1152273d9f13SBarry Smith if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1153da3a660dSBarry Smith for (i=0; i<n; i++) { 1154da3a660dSBarry Smith x = r[i]; 1155da3a660dSBarry Smith v = mat->v + i*m; 1156da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1157da3a660dSBarry Smith } 11581ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1159*efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1160da3a660dSBarry Smith } 11613a40ed3dSBarry Smith PetscFunctionReturn(0); 1162289bc588SBarry Smith } 1163289bc588SBarry Smith 11644a2ae208SSatish Balay #undef __FUNCT__ 11654a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1166dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1167289bc588SBarry Smith { 1168c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 116987828ca2SBarry Smith PetscScalar *v = mat->v; 1170329f5518SBarry Smith PetscReal sum = 0.0; 117113f74950SBarry Smith PetscInt lda=mat->lda,m=A->m,i,j; 1172*efee365bSSatish Balay PetscErrorCode ierr; 117355659b69SBarry Smith 11743a40ed3dSBarry Smith PetscFunctionBegin; 1175289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1176a5ce6ee0Svictorle if (lda>m) { 1177a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1178a5ce6ee0Svictorle v = mat->v+j*lda; 1179a5ce6ee0Svictorle for (i=0; i<m; i++) { 1180a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1181a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1182a5ce6ee0Svictorle #else 1183a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1184a5ce6ee0Svictorle #endif 1185a5ce6ee0Svictorle } 1186a5ce6ee0Svictorle } 1187a5ce6ee0Svictorle } else { 1188273d9f13SBarry Smith for (i=0; i<A->n*A->m; i++) { 1189aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1190329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1191289bc588SBarry Smith #else 1192289bc588SBarry Smith sum += (*v)*(*v); v++; 1193289bc588SBarry Smith #endif 1194289bc588SBarry Smith } 1195a5ce6ee0Svictorle } 1196064f8208SBarry Smith *nrm = sqrt(sum); 1197*efee365bSSatish Balay ierr = PetscLogFlops(2*A->n*A->m);CHKERRQ(ierr); 11983a40ed3dSBarry Smith } else if (type == NORM_1) { 1199064f8208SBarry Smith *nrm = 0.0; 1200273d9f13SBarry Smith for (j=0; j<A->n; j++) { 12011b807ce4Svictorle v = mat->v + j*mat->lda; 1202289bc588SBarry Smith sum = 0.0; 1203273d9f13SBarry Smith for (i=0; i<A->m; i++) { 120433a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1205289bc588SBarry Smith } 1206064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1207289bc588SBarry Smith } 1208*efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 12093a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1210064f8208SBarry Smith *nrm = 0.0; 1211273d9f13SBarry Smith for (j=0; j<A->m; j++) { 1212289bc588SBarry Smith v = mat->v + j; 1213289bc588SBarry Smith sum = 0.0; 1214273d9f13SBarry Smith for (i=0; i<A->n; i++) { 12151b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1216289bc588SBarry Smith } 1217064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1218289bc588SBarry Smith } 1219*efee365bSSatish Balay ierr = PetscLogFlops(A->n*A->m);CHKERRQ(ierr); 12203a40ed3dSBarry Smith } else { 122129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1222289bc588SBarry Smith } 12233a40ed3dSBarry Smith PetscFunctionReturn(0); 1224289bc588SBarry Smith } 1225289bc588SBarry Smith 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1228dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op) 1229289bc588SBarry Smith { 1230c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 123167e560aaSBarry Smith 12323a40ed3dSBarry Smith PetscFunctionBegin; 1233b5a2b587SKris Buschelman switch (op) { 1234b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 1235b5a2b587SKris Buschelman aij->roworiented = PETSC_TRUE; 1236b5a2b587SKris Buschelman break; 1237b5a2b587SKris Buschelman case MAT_COLUMN_ORIENTED: 1238b5a2b587SKris Buschelman aij->roworiented = PETSC_FALSE; 1239b5a2b587SKris Buschelman break; 1240b5a2b587SKris Buschelman case MAT_ROWS_SORTED: 1241b5a2b587SKris Buschelman case MAT_ROWS_UNSORTED: 1242b5a2b587SKris Buschelman case MAT_COLUMNS_SORTED: 1243b5a2b587SKris Buschelman case MAT_COLUMNS_UNSORTED: 1244b5a2b587SKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1245b5a2b587SKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1246b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1247b5a2b587SKris Buschelman case MAT_NO_NEW_DIAGONALS: 1248b5a2b587SKris Buschelman case MAT_YES_NEW_DIAGONALS: 1249b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1250b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 1251b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1252b5a2b587SKris Buschelman break; 125377e54ba9SKris Buschelman case MAT_SYMMETRIC: 125477e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 12559a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 12569a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 12579a4540c5SBarry Smith case MAT_HERMITIAN: 12589a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 12599a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 12609a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 126177e54ba9SKris Buschelman break; 1262b5a2b587SKris Buschelman default: 126329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 12643a40ed3dSBarry Smith } 12653a40ed3dSBarry Smith PetscFunctionReturn(0); 1266289bc588SBarry Smith } 1267289bc588SBarry Smith 12684a2ae208SSatish Balay #undef __FUNCT__ 12694a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1270dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 12716f0a148fSBarry Smith { 1272ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 12736849ba73SBarry Smith PetscErrorCode ierr; 127413f74950SBarry Smith PetscInt lda=l->lda,m=A->m,j; 12753a40ed3dSBarry Smith 12763a40ed3dSBarry Smith PetscFunctionBegin; 1277a5ce6ee0Svictorle if (lda>m) { 1278a5ce6ee0Svictorle for (j=0; j<A->n; j++) { 1279a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1280a5ce6ee0Svictorle } 1281a5ce6ee0Svictorle } else { 128287828ca2SBarry Smith ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1283a5ce6ee0Svictorle } 12843a40ed3dSBarry Smith PetscFunctionReturn(0); 12856f0a148fSBarry Smith } 12866f0a148fSBarry Smith 12874a2ae208SSatish Balay #undef __FUNCT__ 12884a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1289dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag) 12906f0a148fSBarry Smith { 1291ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 12926849ba73SBarry Smith PetscErrorCode ierr; 129313f74950SBarry Smith PetscInt n = A->n,i,j,N,*rows; 129487828ca2SBarry Smith PetscScalar *slot; 129555659b69SBarry Smith 12963a40ed3dSBarry Smith PetscFunctionBegin; 1297e03a110bSBarry Smith ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 129878b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 12996f0a148fSBarry Smith for (i=0; i<N; i++) { 13006f0a148fSBarry Smith slot = l->v + rows[i]; 13016f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13026f0a148fSBarry Smith } 13036f0a148fSBarry Smith if (diag) { 13046f0a148fSBarry Smith for (i=0; i<N; i++) { 13056f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1306260925b8SBarry Smith *slot = *diag; 13076f0a148fSBarry Smith } 13086f0a148fSBarry Smith } 1309260925b8SBarry Smith ISRestoreIndices(is,&rows); 13103a40ed3dSBarry Smith PetscFunctionReturn(0); 13116f0a148fSBarry Smith } 1312557bce09SLois Curfman McInnes 13134a2ae208SSatish Balay #undef __FUNCT__ 13144a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1315dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 131664e87e97SBarry Smith { 1317c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13183a40ed3dSBarry Smith 13193a40ed3dSBarry Smith PetscFunctionBegin; 132064e87e97SBarry Smith *array = mat->v; 13213a40ed3dSBarry Smith PetscFunctionReturn(0); 132264e87e97SBarry Smith } 13230754003eSLois Curfman McInnes 13244a2ae208SSatish Balay #undef __FUNCT__ 13254a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1326dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1327ff14e315SSatish Balay { 13283a40ed3dSBarry Smith PetscFunctionBegin; 132909b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13303a40ed3dSBarry Smith PetscFunctionReturn(0); 1331ff14e315SSatish Balay } 13320754003eSLois Curfman McInnes 13334a2ae208SSatish Balay #undef __FUNCT__ 13344a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 133513f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 13360754003eSLois Curfman McInnes { 1337c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13386849ba73SBarry Smith PetscErrorCode ierr; 133913f74950SBarry Smith PetscInt i,j,m = A->m,*irow,*icol,nrows,ncols; 134087828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13410754003eSLois Curfman McInnes Mat newmat; 13420754003eSLois Curfman McInnes 13433a40ed3dSBarry Smith PetscFunctionBegin; 134478b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 134578b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1346e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1347e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13480754003eSLois Curfman McInnes 1349182d2002SSatish Balay /* Check submatrixcall */ 1350182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 135113f74950SBarry Smith PetscInt n_cols,n_rows; 1352182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 135329bbc08cSBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1354182d2002SSatish Balay newmat = *B; 1355182d2002SSatish Balay } else { 13560754003eSLois Curfman McInnes /* Create and fill new matrix */ 13575c5985e7SKris Buschelman ierr = MatCreate(A->comm,nrows,ncols,nrows,ncols,&newmat);CHKERRQ(ierr); 13585c5985e7SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 13595c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1360182d2002SSatish Balay } 1361182d2002SSatish Balay 1362182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1363182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1364182d2002SSatish Balay 1365182d2002SSatish Balay for (i=0; i<ncols; i++) { 1366182d2002SSatish Balay av = v + m*icol[i]; 1367182d2002SSatish Balay for (j=0; j<nrows; j++) { 1368182d2002SSatish Balay *bv++ = av[irow[j]]; 13690754003eSLois Curfman McInnes } 13700754003eSLois Curfman McInnes } 1371182d2002SSatish Balay 1372182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13736d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13746d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13750754003eSLois Curfman McInnes 13760754003eSLois Curfman McInnes /* Free work space */ 137778b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 137878b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1379182d2002SSatish Balay *B = newmat; 13803a40ed3dSBarry Smith PetscFunctionReturn(0); 13810754003eSLois Curfman McInnes } 13820754003eSLois Curfman McInnes 13834a2ae208SSatish Balay #undef __FUNCT__ 13844a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 138513f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1386905e6a2fSBarry Smith { 13876849ba73SBarry Smith PetscErrorCode ierr; 138813f74950SBarry Smith PetscInt i; 1389905e6a2fSBarry Smith 13903a40ed3dSBarry Smith PetscFunctionBegin; 1391905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1392b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1393905e6a2fSBarry Smith } 1394905e6a2fSBarry Smith 1395905e6a2fSBarry Smith for (i=0; i<n; i++) { 13966a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1397905e6a2fSBarry Smith } 13983a40ed3dSBarry Smith PetscFunctionReturn(0); 1399905e6a2fSBarry Smith } 1400905e6a2fSBarry Smith 14014a2ae208SSatish Balay #undef __FUNCT__ 14024a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1403dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14044b0e389bSBarry Smith { 14054b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 14066849ba73SBarry Smith PetscErrorCode ierr; 140713f74950SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j; 14083a40ed3dSBarry Smith 14093a40ed3dSBarry Smith PetscFunctionBegin; 141033f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 141133f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1412cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 14133a40ed3dSBarry Smith PetscFunctionReturn(0); 14143a40ed3dSBarry Smith } 14150dbb7854Svictorle if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1416a5ce6ee0Svictorle if (lda1>m || lda2>m) { 14170dbb7854Svictorle for (j=0; j<n; j++) { 1418a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1419a5ce6ee0Svictorle } 1420a5ce6ee0Svictorle } else { 142187828ca2SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1422a5ce6ee0Svictorle } 1423273d9f13SBarry Smith PetscFunctionReturn(0); 1424273d9f13SBarry Smith } 1425273d9f13SBarry Smith 14264a2ae208SSatish Balay #undef __FUNCT__ 14274a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1428dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1429273d9f13SBarry Smith { 1430dfbe8321SBarry Smith PetscErrorCode ierr; 1431273d9f13SBarry Smith 1432273d9f13SBarry Smith PetscFunctionBegin; 1433273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 14343a40ed3dSBarry Smith PetscFunctionReturn(0); 14354b0e389bSBarry Smith } 14364b0e389bSBarry Smith 1437289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1438a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1439905e6a2fSBarry Smith MatGetRow_SeqDense, 1440905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1441905e6a2fSBarry Smith MatMult_SeqDense, 144297304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 14437c922b88SBarry Smith MatMultTranspose_SeqDense, 14447c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1445905e6a2fSBarry Smith MatSolve_SeqDense, 1446905e6a2fSBarry Smith MatSolveAdd_SeqDense, 14477c922b88SBarry Smith MatSolveTranspose_SeqDense, 144897304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense, 1449905e6a2fSBarry Smith MatLUFactor_SeqDense, 1450905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1451ec8511deSBarry Smith MatRelax_SeqDense, 1452ec8511deSBarry Smith MatTranspose_SeqDense, 145397304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1454905e6a2fSBarry Smith MatEqual_SeqDense, 1455905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1456905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1457905e6a2fSBarry Smith MatNorm_SeqDense, 145897304618SKris Buschelman /*20*/ 0, 1459905e6a2fSBarry Smith 0, 1460905e6a2fSBarry Smith 0, 1461905e6a2fSBarry Smith MatSetOption_SeqDense, 1462905e6a2fSBarry Smith MatZeroEntries_SeqDense, 146397304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1464905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1465905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1466905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1467905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 146897304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1469273d9f13SBarry Smith 0, 1470905e6a2fSBarry Smith 0, 1471905e6a2fSBarry Smith MatGetArray_SeqDense, 1472905e6a2fSBarry Smith MatRestoreArray_SeqDense, 147397304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1474a5ae1ecdSBarry Smith 0, 1475a5ae1ecdSBarry Smith 0, 1476a5ae1ecdSBarry Smith 0, 1477a5ae1ecdSBarry Smith 0, 147897304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1479a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1480a5ae1ecdSBarry Smith 0, 14814b0e389bSBarry Smith MatGetValues_SeqDense, 1482a5ae1ecdSBarry Smith MatCopy_SeqDense, 148397304618SKris Buschelman /*45*/ 0, 1484a5ae1ecdSBarry Smith MatScale_SeqDense, 1485a5ae1ecdSBarry Smith 0, 1486a5ae1ecdSBarry Smith 0, 1487a5ae1ecdSBarry Smith 0, 1488521d7252SBarry Smith /*50*/ 0, 1489a5ae1ecdSBarry Smith 0, 1490a5ae1ecdSBarry Smith 0, 1491a5ae1ecdSBarry Smith 0, 1492a5ae1ecdSBarry Smith 0, 149397304618SKris Buschelman /*55*/ 0, 1494a5ae1ecdSBarry Smith 0, 1495a5ae1ecdSBarry Smith 0, 1496a5ae1ecdSBarry Smith 0, 1497a5ae1ecdSBarry Smith 0, 149897304618SKris Buschelman /*60*/ 0, 1499e03a110bSBarry Smith MatDestroy_SeqDense, 1500e03a110bSBarry Smith MatView_SeqDense, 150197304618SKris Buschelman MatGetPetscMaps_Petsc, 150297304618SKris Buschelman 0, 150397304618SKris Buschelman /*65*/ 0, 150497304618SKris Buschelman 0, 150597304618SKris Buschelman 0, 150697304618SKris Buschelman 0, 150797304618SKris Buschelman 0, 150897304618SKris Buschelman /*70*/ 0, 150997304618SKris Buschelman 0, 151097304618SKris Buschelman 0, 151197304618SKris Buschelman 0, 151297304618SKris Buschelman 0, 151397304618SKris Buschelman /*75*/ 0, 151497304618SKris Buschelman 0, 151597304618SKris Buschelman 0, 151697304618SKris Buschelman 0, 151797304618SKris Buschelman 0, 151897304618SKris Buschelman /*80*/ 0, 151997304618SKris Buschelman 0, 152097304618SKris Buschelman 0, 152197304618SKris Buschelman 0, 1522865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense, 1523865e5f61SKris Buschelman 0, 1524865e5f61SKris Buschelman 0, 1525865e5f61SKris Buschelman 0, 1526865e5f61SKris Buschelman 0, 1527865e5f61SKris Buschelman 0, 1528865e5f61SKris Buschelman /*90*/ 0, 1529865e5f61SKris Buschelman 0, 1530865e5f61SKris Buschelman 0, 1531865e5f61SKris Buschelman 0, 1532865e5f61SKris Buschelman 0, 1533865e5f61SKris Buschelman /*95*/ 0, 1534865e5f61SKris Buschelman 0, 1535865e5f61SKris Buschelman 0, 1536865e5f61SKris Buschelman 0}; 153790ace30eSBarry Smith 15384a2ae208SSatish Balay #undef __FUNCT__ 15394a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 15404b828684SBarry Smith /*@C 1541fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1542d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1543d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1544289bc588SBarry Smith 1545db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1546db81eaa0SLois Curfman McInnes 154720563c6bSBarry Smith Input Parameters: 1548db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 15490c775827SLois Curfman McInnes . m - number of rows 155018f449edSLois Curfman McInnes . n - number of columns 1551db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1552dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 155320563c6bSBarry Smith 155420563c6bSBarry Smith Output Parameter: 155544cd7ae7SLois Curfman McInnes . A - the matrix 155620563c6bSBarry Smith 1557b259b22eSLois Curfman McInnes Notes: 155818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 155918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1560b4fd4287SBarry Smith set data=PETSC_NULL. 156118f449edSLois Curfman McInnes 1562027ccd11SLois Curfman McInnes Level: intermediate 1563027ccd11SLois Curfman McInnes 1564dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1565d65003e9SLois Curfman McInnes 1566db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 156720563c6bSBarry Smith @*/ 156813f74950SBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1569289bc588SBarry Smith { 1570dfbe8321SBarry Smith PetscErrorCode ierr; 15713b2fbd54SBarry Smith 15723a40ed3dSBarry Smith PetscFunctionBegin; 1573273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1574273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1575273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1576273d9f13SBarry Smith PetscFunctionReturn(0); 1577273d9f13SBarry Smith } 1578273d9f13SBarry Smith 15794a2ae208SSatish Balay #undef __FUNCT__ 15804a2ae208SSatish Balay #define __FUNCT__ "MatSeqDensePreallocation" 1581273d9f13SBarry Smith /*@C 1582273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1583273d9f13SBarry Smith 1584273d9f13SBarry Smith Collective on MPI_Comm 1585273d9f13SBarry Smith 1586273d9f13SBarry Smith Input Parameters: 1587273d9f13SBarry Smith + A - the matrix 1588273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1589273d9f13SBarry Smith 1590273d9f13SBarry Smith Notes: 1591273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1592273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1593273d9f13SBarry Smith set data=PETSC_NULL. 1594273d9f13SBarry Smith 1595273d9f13SBarry Smith Level: intermediate 1596273d9f13SBarry Smith 1597273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1598273d9f13SBarry Smith 1599273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1600273d9f13SBarry Smith @*/ 1601dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1602273d9f13SBarry Smith { 1603dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1604a23d5eceSKris Buschelman 1605a23d5eceSKris Buschelman PetscFunctionBegin; 1606a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1607a23d5eceSKris Buschelman if (f) { 1608a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1609a23d5eceSKris Buschelman } 1610a23d5eceSKris Buschelman PetscFunctionReturn(0); 1611a23d5eceSKris Buschelman } 1612a23d5eceSKris Buschelman 1613a23d5eceSKris Buschelman EXTERN_C_BEGIN 1614a23d5eceSKris Buschelman #undef __FUNCT__ 1615a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1616dfbe8321SBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1617a23d5eceSKris Buschelman { 1618273d9f13SBarry Smith Mat_SeqDense *b; 1619dfbe8321SBarry Smith PetscErrorCode ierr; 1620273d9f13SBarry Smith 1621273d9f13SBarry Smith PetscFunctionBegin; 1622273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1623273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1624273d9f13SBarry Smith if (!data) { 162587828ca2SBarry Smith ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 162687828ca2SBarry Smith ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1627273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 162852e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar));CHKERRQ(ierr); 1629273d9f13SBarry Smith } else { /* user-allocated storage */ 1630273d9f13SBarry Smith b->v = data; 1631273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1632273d9f13SBarry Smith } 1633273d9f13SBarry Smith PetscFunctionReturn(0); 1634273d9f13SBarry Smith } 1635a23d5eceSKris Buschelman EXTERN_C_END 1636273d9f13SBarry Smith 16371b807ce4Svictorle #undef __FUNCT__ 16381b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 16391b807ce4Svictorle /*@C 16401b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 16411b807ce4Svictorle 16421b807ce4Svictorle Input parameter: 16431b807ce4Svictorle + A - the matrix 16441b807ce4Svictorle - lda - the leading dimension 16451b807ce4Svictorle 16461b807ce4Svictorle Notes: 16471b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 16481b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 16491b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 16501b807ce4Svictorle 16511b807ce4Svictorle Level: intermediate 16521b807ce4Svictorle 16531b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 16541b807ce4Svictorle 16551b807ce4Svictorle .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 16561b807ce4Svictorle @*/ 165713f74950SBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 16581b807ce4Svictorle { 16591b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 16601b807ce4Svictorle PetscFunctionBegin; 166177431f27SBarry Smith if (lda < B->m) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->m); 16621b807ce4Svictorle b->lda = lda; 16631b807ce4Svictorle PetscFunctionReturn(0); 16641b807ce4Svictorle } 16651b807ce4Svictorle 16660bad9183SKris Buschelman /*MC 1667fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 16680bad9183SKris Buschelman 16690bad9183SKris Buschelman Options Database Keys: 16700bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 16710bad9183SKris Buschelman 16720bad9183SKris Buschelman Level: beginner 16730bad9183SKris Buschelman 16740bad9183SKris Buschelman .seealso: MatCreateSeqDense 16750bad9183SKris Buschelman M*/ 16760bad9183SKris Buschelman 1677273d9f13SBarry Smith EXTERN_C_BEGIN 16784a2ae208SSatish Balay #undef __FUNCT__ 16794a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1680dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 1681273d9f13SBarry Smith { 1682273d9f13SBarry Smith Mat_SeqDense *b; 1683dfbe8321SBarry Smith PetscErrorCode ierr; 16847c334f02SBarry Smith PetscMPIInt size; 1685273d9f13SBarry Smith 1686273d9f13SBarry Smith PetscFunctionBegin; 1687273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 168829bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 168955659b69SBarry Smith 1690273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1691273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1692273d9f13SBarry Smith 1693b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1694549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 169544cd7ae7SLois Curfman McInnes B->factor = 0; 169690f02eecSBarry Smith B->mapping = 0; 169752e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr); 169844cd7ae7SLois Curfman McInnes B->data = (void*)b; 169918f449edSLois Curfman McInnes 17008a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 17018a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1702a5ae1ecdSBarry Smith 170344cd7ae7SLois Curfman McInnes b->pivots = 0; 1704273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1705273d9f13SBarry Smith b->v = 0; 17061b807ce4Svictorle b->lda = B->m; 17074e220ebcSLois Curfman McInnes 1708a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1709a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 1710a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 17113a40ed3dSBarry Smith PetscFunctionReturn(0); 1712289bc588SBarry Smith } 1713273d9f13SBarry Smith EXTERN_C_END 1714