1be1d678aSKris Buschelman 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 6dec5eb66SMatthew G Knepley #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <petscblaslapack.h> 8289bc588SBarry Smith 96a63e612SBarry Smith #include <../src/mat/impls/aij/seq/aij.h> 10b2573a8aSBarry Smith 116a63e612SBarry Smith #undef __FUNCT__ 126a63e612SBarry Smith #define __FUNCT__ "MatConvert_SeqDense_SeqAIJ" 138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 146a63e612SBarry Smith { 156a63e612SBarry Smith Mat B; 166a63e612SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 176a63e612SBarry Smith PetscErrorCode ierr; 189399e1b8SMatthew G. Knepley PetscInt i, j; 199399e1b8SMatthew G. Knepley PetscInt *rows, *nnz; 209399e1b8SMatthew G. Knepley MatScalar *aa = a->v, *vals; 216a63e612SBarry Smith 226a63e612SBarry Smith PetscFunctionBegin; 23ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 246a63e612SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 256a63e612SBarry Smith ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 269399e1b8SMatthew G. Knepley ierr = PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals);CHKERRQ(ierr); 279399e1b8SMatthew G. Knepley for (j=0; j<A->cmap->n; j++) { 289399e1b8SMatthew G. Knepley for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) ++nnz[i]; 296a63e612SBarry Smith aa += a->lda; 306a63e612SBarry Smith } 319399e1b8SMatthew G. Knepley ierr = MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz);CHKERRQ(ierr); 329399e1b8SMatthew G. Knepley aa = a->v; 339399e1b8SMatthew G. Knepley for (j=0; j<A->cmap->n; j++) { 349399e1b8SMatthew G. Knepley PetscInt numRows = 0; 359399e1b8SMatthew G. Knepley for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || i == j) {rows[numRows] = i; vals[numRows++] = aa[i];} 369399e1b8SMatthew G. Knepley ierr = MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES);CHKERRQ(ierr); 379399e1b8SMatthew G. Knepley aa += a->lda; 389399e1b8SMatthew G. Knepley } 399399e1b8SMatthew G. Knepley ierr = PetscFree3(rows,nnz,vals);CHKERRQ(ierr); 406a63e612SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 416a63e612SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 426a63e612SBarry Smith 436a63e612SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 446a63e612SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 456a63e612SBarry Smith } else { 466a63e612SBarry Smith *newmat = B; 476a63e612SBarry Smith } 486a63e612SBarry Smith PetscFunctionReturn(0); 496a63e612SBarry Smith } 506a63e612SBarry Smith 514a2ae208SSatish Balay #undef __FUNCT__ 524a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 53f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 541987afe7SBarry Smith { 551987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 56f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 5713f74950SBarry Smith PetscInt j; 580805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 59efee365bSSatish Balay PetscErrorCode ierr; 603a40ed3dSBarry Smith 613a40ed3dSBarry Smith PetscFunctionBegin; 62c5df96a5SBarry Smith ierr = PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);CHKERRQ(ierr); 63c5df96a5SBarry Smith ierr = PetscBLASIntCast(X->rmap->n,&m);CHKERRQ(ierr); 64c5df96a5SBarry Smith ierr = PetscBLASIntCast(x->lda,&ldax);CHKERRQ(ierr); 65c5df96a5SBarry Smith ierr = PetscBLASIntCast(y->lda,&lday);CHKERRQ(ierr); 66a5ce6ee0Svictorle if (ldax>m || lday>m) { 67d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 688b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one)); 69a5ce6ee0Svictorle } 70a5ce6ee0Svictorle } else { 718b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one)); 72a5ce6ee0Svictorle } 730450473dSBarry Smith ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr); 743a40ed3dSBarry Smith PetscFunctionReturn(0); 751987afe7SBarry Smith } 761987afe7SBarry Smith 774a2ae208SSatish Balay #undef __FUNCT__ 784a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 79dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 80289bc588SBarry Smith { 81d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 823a40ed3dSBarry Smith 833a40ed3dSBarry Smith PetscFunctionBegin; 844e220ebcSLois Curfman McInnes info->block_size = 1.0; 854e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 866de62eeeSBarry Smith info->nz_used = (double)N; 876de62eeeSBarry Smith info->nz_unneeded = (double)0; 884e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 894e220ebcSLois Curfman McInnes info->mallocs = 0; 907adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 914e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 924e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 934e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 943a40ed3dSBarry Smith PetscFunctionReturn(0); 95289bc588SBarry Smith } 96289bc588SBarry Smith 974a2ae208SSatish Balay #undef __FUNCT__ 984a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 99f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 10080cd9d93SLois Curfman McInnes { 101273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 102f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 103efee365bSSatish Balay PetscErrorCode ierr; 104c5df96a5SBarry Smith PetscBLASInt one = 1,j,nz,lda; 10580cd9d93SLois Curfman McInnes 1063a40ed3dSBarry Smith PetscFunctionBegin; 107c5df96a5SBarry Smith ierr = PetscBLASIntCast(a->lda,&lda);CHKERRQ(ierr); 108d0f46423SBarry Smith if (lda>A->rmap->n) { 109c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&nz);CHKERRQ(ierr); 110d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1118b83055fSJed Brown PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v+j*lda,&one)); 112a5ce6ee0Svictorle } 113a5ce6ee0Svictorle } else { 114c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);CHKERRQ(ierr); 1158b83055fSJed Brown PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&oalpha,a->v,&one)); 116a5ce6ee0Svictorle } 117efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 1183a40ed3dSBarry Smith PetscFunctionReturn(0); 11980cd9d93SLois Curfman McInnes } 12080cd9d93SLois Curfman McInnes 1211cbb95d3SBarry Smith #undef __FUNCT__ 1221cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 123ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 1241cbb95d3SBarry Smith { 1251cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 126d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,N; 1271cbb95d3SBarry Smith PetscScalar *v = a->v; 1281cbb95d3SBarry Smith 1291cbb95d3SBarry Smith PetscFunctionBegin; 1301cbb95d3SBarry Smith *fl = PETSC_FALSE; 131d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 1321cbb95d3SBarry Smith N = a->lda; 1331cbb95d3SBarry Smith 1341cbb95d3SBarry Smith for (i=0; i<m; i++) { 1351cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 1361cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 1371cbb95d3SBarry Smith } 1381cbb95d3SBarry Smith } 1391cbb95d3SBarry Smith *fl = PETSC_TRUE; 1401cbb95d3SBarry Smith PetscFunctionReturn(0); 1411cbb95d3SBarry Smith } 1421cbb95d3SBarry Smith 143b24902e0SBarry Smith #undef __FUNCT__ 144b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 145719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 146b24902e0SBarry Smith { 147b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 148b24902e0SBarry Smith PetscErrorCode ierr; 149b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 150b24902e0SBarry Smith 151b24902e0SBarry Smith PetscFunctionBegin; 152aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr); 153aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr); 1540298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,NULL);CHKERRQ(ierr); 155b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 156b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 157d0f46423SBarry Smith if (lda>A->rmap->n) { 158d0f46423SBarry Smith m = A->rmap->n; 159d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 160b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 161b24902e0SBarry Smith } 162b24902e0SBarry Smith } else { 163d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 164b24902e0SBarry Smith } 165b24902e0SBarry Smith } 166b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 167b24902e0SBarry Smith PetscFunctionReturn(0); 168b24902e0SBarry Smith } 169b24902e0SBarry Smith 1704a2ae208SSatish Balay #undef __FUNCT__ 1714a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 172dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 17302cad45dSBarry Smith { 1746849ba73SBarry Smith PetscErrorCode ierr; 17502cad45dSBarry Smith 1763a40ed3dSBarry Smith PetscFunctionBegin; 177ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),newmat);CHKERRQ(ierr); 178d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1795c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 180719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 181b24902e0SBarry Smith PetscFunctionReturn(0); 182b24902e0SBarry Smith } 183b24902e0SBarry Smith 1846ee01492SSatish Balay 1850481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 186719d5645SBarry Smith 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1890481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 190289bc588SBarry Smith { 1914482741eSBarry Smith MatFactorInfo info; 192a093e273SMatthew Knepley PetscErrorCode ierr; 1933a40ed3dSBarry Smith 1943a40ed3dSBarry Smith PetscFunctionBegin; 195c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 196719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1973a40ed3dSBarry Smith PetscFunctionReturn(0); 198289bc588SBarry Smith } 1996ee01492SSatish Balay 2000b4b3355SBarry Smith #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 202dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 203289bc588SBarry Smith { 204c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2056849ba73SBarry Smith PetscErrorCode ierr; 20687828ca2SBarry Smith PetscScalar *x,*y; 207c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 20867e560aaSBarry Smith 2093a40ed3dSBarry Smith PetscFunctionBegin; 210c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2111ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2121ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 213d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 214d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 215ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 216e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 217ae7cfcebSSatish Balay #else 2188b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 219e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 220ae7cfcebSSatish Balay #endif 221d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_CHOLESKY) { 222ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 223e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 224ae7cfcebSSatish Balay #else 2258b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 226e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 227ae7cfcebSSatish Balay #endif 2282205254eSKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2291ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2301ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 231dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2323a40ed3dSBarry Smith PetscFunctionReturn(0); 233289bc588SBarry Smith } 2346ee01492SSatish Balay 2354a2ae208SSatish Balay #undef __FUNCT__ 23685e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense" 23785e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X) 23885e2c93fSHong Zhang { 23985e2c93fSHong Zhang Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 24085e2c93fSHong Zhang PetscErrorCode ierr; 24185e2c93fSHong Zhang PetscScalar *b,*x; 242efb80c78SLisandro Dalcin PetscInt n; 243783b601eSJed Brown PetscBLASInt nrhs,info,m; 244bda8bf91SBarry Smith PetscBool flg; 24585e2c93fSHong Zhang 24685e2c93fSHong Zhang PetscFunctionBegin; 247c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2480298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 249ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 2500298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 251ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 252bda8bf91SBarry Smith 2530298fd71SBarry Smith ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr); 254c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr); 2558c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 2568c778c55SBarry Smith ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr); 25785e2c93fSHong Zhang 258f272c775SHong Zhang ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 25985e2c93fSHong Zhang 26085e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 26185e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 26285e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 26385e2c93fSHong Zhang #else 2648b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 26585e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 26685e2c93fSHong Zhang #endif 26785e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY) { 26885e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 26985e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 27085e2c93fSHong Zhang #else 2718b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info)); 27285e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 27385e2c93fSHong Zhang #endif 2742205254eSKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 27585e2c93fSHong Zhang 2768c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2778c778c55SBarry Smith ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr); 27885e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr); 27985e2c93fSHong Zhang PetscFunctionReturn(0); 28085e2c93fSHong Zhang } 28185e2c93fSHong Zhang 28285e2c93fSHong Zhang #undef __FUNCT__ 2834a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 284dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 285da3a660dSBarry Smith { 286c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 287dfbe8321SBarry Smith PetscErrorCode ierr; 28887828ca2SBarry Smith PetscScalar *x,*y; 289c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 29067e560aaSBarry Smith 2913a40ed3dSBarry Smith PetscFunctionBegin; 292c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2931ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2941ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 295d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 296752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 297da3a660dSBarry Smith if (mat->pivots) { 298ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 299e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 300ae7cfcebSSatish Balay #else 3018b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 302e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 303ae7cfcebSSatish Balay #endif 3047a97a34bSBarry Smith } else { 305ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 306e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 307ae7cfcebSSatish Balay #else 3088b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 309e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 310ae7cfcebSSatish Balay #endif 311da3a660dSBarry Smith } 3121ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3131ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 314dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 316da3a660dSBarry Smith } 3176ee01492SSatish Balay 3184a2ae208SSatish Balay #undef __FUNCT__ 3194a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 320dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 321da3a660dSBarry Smith { 322c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 323dfbe8321SBarry Smith PetscErrorCode ierr; 32487828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 325da3a660dSBarry Smith Vec tmp = 0; 326c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 32767e560aaSBarry Smith 3283a40ed3dSBarry Smith PetscFunctionBegin; 329c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 3301ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3311ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 332d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 333da3a660dSBarry Smith if (yy == zz) { 33478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 3353bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);CHKERRQ(ierr); 33678b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 337da3a660dSBarry Smith } 338d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 339752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 340da3a660dSBarry Smith if (mat->pivots) { 341ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 342e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 343ae7cfcebSSatish Balay #else 3448b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 345e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 346ae7cfcebSSatish Balay #endif 347a8c6a408SBarry Smith } else { 348ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 349e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 350ae7cfcebSSatish Balay #else 3518b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 352e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 353ae7cfcebSSatish Balay #endif 354da3a660dSBarry Smith } 3556bf464f9SBarry Smith if (tmp) { 3566bf464f9SBarry Smith ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3576bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3586bf464f9SBarry Smith } else { 3596bf464f9SBarry Smith ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 3606bf464f9SBarry Smith } 3611ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3621ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 363dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3643a40ed3dSBarry Smith PetscFunctionReturn(0); 365da3a660dSBarry Smith } 36667e560aaSBarry Smith 3674a2ae208SSatish Balay #undef __FUNCT__ 3684a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 369dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 370da3a660dSBarry Smith { 371c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3726849ba73SBarry Smith PetscErrorCode ierr; 37387828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 374da3a660dSBarry Smith Vec tmp; 375c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 37667e560aaSBarry Smith 3773a40ed3dSBarry Smith PetscFunctionBegin; 378c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 379d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3801ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3811ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 382da3a660dSBarry Smith if (yy == zz) { 38378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 3843bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)tmp);CHKERRQ(ierr); 38578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 386da3a660dSBarry Smith } 387d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 388752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 389da3a660dSBarry Smith if (mat->pivots) { 390ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 391e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 392ae7cfcebSSatish Balay #else 3938b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 394e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 395ae7cfcebSSatish Balay #endif 3963a40ed3dSBarry Smith } else { 397ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 398e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 399ae7cfcebSSatish Balay #else 4008b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 401e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 402ae7cfcebSSatish Balay #endif 403da3a660dSBarry Smith } 40490f02eecSBarry Smith if (tmp) { 4052dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 4066bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 4073a40ed3dSBarry Smith } else { 4082dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 40990f02eecSBarry Smith } 4101ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4111ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 412dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 4133a40ed3dSBarry Smith PetscFunctionReturn(0); 414da3a660dSBarry Smith } 415db4efbfdSBarry Smith 416db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 417db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 418db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 419db4efbfdSBarry Smith #undef __FUNCT__ 420db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 4210481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 422db4efbfdSBarry Smith { 423db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 424db4efbfdSBarry Smith PetscFunctionBegin; 425e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 426db4efbfdSBarry Smith #else 427db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 428db4efbfdSBarry Smith PetscErrorCode ierr; 429db4efbfdSBarry Smith PetscBLASInt n,m,info; 430db4efbfdSBarry Smith 431db4efbfdSBarry Smith PetscFunctionBegin; 432c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 433c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 434db4efbfdSBarry Smith if (!mat->pivots) { 435785e854fSJed Brown ierr = PetscMalloc1((A->rmap->n+1),&mat->pivots);CHKERRQ(ierr); 4363bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 437db4efbfdSBarry Smith } 438db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 4398e57ea43SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 4408b83055fSJed Brown PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info)); 4418e57ea43SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 4428e57ea43SSatish Balay 443e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 444e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 445db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 446db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 447db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 448db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 449d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 450db4efbfdSBarry Smith 451dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 452db4efbfdSBarry Smith #endif 453db4efbfdSBarry Smith PetscFunctionReturn(0); 454db4efbfdSBarry Smith } 455db4efbfdSBarry Smith 456db4efbfdSBarry Smith #undef __FUNCT__ 457db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4580481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 459db4efbfdSBarry Smith { 460db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 461db4efbfdSBarry Smith PetscFunctionBegin; 462e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 463db4efbfdSBarry Smith #else 464db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 465db4efbfdSBarry Smith PetscErrorCode ierr; 466c5df96a5SBarry Smith PetscBLASInt info,n; 467db4efbfdSBarry Smith 468db4efbfdSBarry Smith PetscFunctionBegin; 469c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 470db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 471db4efbfdSBarry Smith 472db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 4738b83055fSJed Brown PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info)); 474e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 475db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 476db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 477db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 478db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 479d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 4802205254eSKarl Rupp 481dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 482db4efbfdSBarry Smith #endif 483db4efbfdSBarry Smith PetscFunctionReturn(0); 484db4efbfdSBarry Smith } 485db4efbfdSBarry Smith 486db4efbfdSBarry Smith 487db4efbfdSBarry Smith #undef __FUNCT__ 488db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4890481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 490db4efbfdSBarry Smith { 491db4efbfdSBarry Smith PetscErrorCode ierr; 492db4efbfdSBarry Smith MatFactorInfo info; 493db4efbfdSBarry Smith 494db4efbfdSBarry Smith PetscFunctionBegin; 495db4efbfdSBarry Smith info.fill = 1.0; 4962205254eSKarl Rupp 497c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 498719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 499db4efbfdSBarry Smith PetscFunctionReturn(0); 500db4efbfdSBarry Smith } 501db4efbfdSBarry Smith 502db4efbfdSBarry Smith #undef __FUNCT__ 503db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 5040481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 505db4efbfdSBarry Smith { 506db4efbfdSBarry Smith PetscFunctionBegin; 507c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 5081bbcc794SSatish Balay fact->preallocated = PETSC_TRUE; 509719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 510db4efbfdSBarry Smith PetscFunctionReturn(0); 511db4efbfdSBarry Smith } 512db4efbfdSBarry Smith 513db4efbfdSBarry Smith #undef __FUNCT__ 514db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 5150481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 516db4efbfdSBarry Smith { 517db4efbfdSBarry Smith PetscFunctionBegin; 518b66fe19dSMatthew G Knepley fact->preallocated = PETSC_TRUE; 519c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 520719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 521db4efbfdSBarry Smith PetscFunctionReturn(0); 522db4efbfdSBarry Smith } 523db4efbfdSBarry Smith 524db4efbfdSBarry Smith #undef __FUNCT__ 525db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 5268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 527db4efbfdSBarry Smith { 528db4efbfdSBarry Smith PetscErrorCode ierr; 529db4efbfdSBarry Smith 530db4efbfdSBarry Smith PetscFunctionBegin; 531ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),fact);CHKERRQ(ierr); 532db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 533db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 534db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU) { 535db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 536db4efbfdSBarry Smith } else { 537db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 538db4efbfdSBarry Smith } 539d5f3da31SBarry Smith (*fact)->factortype = ftype; 540db4efbfdSBarry Smith PetscFunctionReturn(0); 541db4efbfdSBarry Smith } 542db4efbfdSBarry Smith 543289bc588SBarry Smith /* ------------------------------------------------------------------*/ 5444a2ae208SSatish Balay #undef __FUNCT__ 54541f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 54641f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 547289bc588SBarry Smith { 548c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54987828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 550dfbe8321SBarry Smith PetscErrorCode ierr; 551d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 552c5df96a5SBarry Smith PetscBLASInt o = 1,bm; 553289bc588SBarry Smith 5543a40ed3dSBarry Smith PetscFunctionBegin; 555c5df96a5SBarry Smith ierr = PetscBLASIntCast(m,&bm);CHKERRQ(ierr); 556289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 55771044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 5582dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 559289bc588SBarry Smith } 5601ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5611ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 562b965ef7fSBarry Smith its = its*lits; 563e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 564289bc588SBarry Smith while (its--) { 565fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 566289bc588SBarry Smith for (i=0; i<m; i++) { 5678b83055fSJed Brown PetscStackCallBLAS("BLASdot",xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o)); 56855a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 569289bc588SBarry Smith } 570289bc588SBarry Smith } 571fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 572289bc588SBarry Smith for (i=m-1; i>=0; i--) { 5738b83055fSJed Brown PetscStackCallBLAS("BLASdot",xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o)); 57455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 575289bc588SBarry Smith } 576289bc588SBarry Smith } 577289bc588SBarry Smith } 5781ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5791ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5803a40ed3dSBarry Smith PetscFunctionReturn(0); 581289bc588SBarry Smith } 582289bc588SBarry Smith 583289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5844a2ae208SSatish Balay #undef __FUNCT__ 5854a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 586dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 587289bc588SBarry Smith { 588c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 58987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 590dfbe8321SBarry Smith PetscErrorCode ierr; 5910805154bSBarry Smith PetscBLASInt m, n,_One=1; 592ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5933a40ed3dSBarry Smith 5943a40ed3dSBarry Smith PetscFunctionBegin; 595c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 596c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 597d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5981ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5991ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6008b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One)); 6011ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6021ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 603dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 6043a40ed3dSBarry Smith PetscFunctionReturn(0); 605289bc588SBarry Smith } 606800995b7SMatthew Knepley 6074a2ae208SSatish Balay #undef __FUNCT__ 6084a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 609dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 610289bc588SBarry Smith { 611c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 61287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 613dfbe8321SBarry Smith PetscErrorCode ierr; 6140805154bSBarry Smith PetscBLASInt m, n, _One=1; 6153a40ed3dSBarry Smith 6163a40ed3dSBarry Smith PetscFunctionBegin; 617c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 618c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 619d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 6201ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6211ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6228b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One)); 6231ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6241ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 625dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 6263a40ed3dSBarry Smith PetscFunctionReturn(0); 627289bc588SBarry Smith } 6286ee01492SSatish Balay 6294a2ae208SSatish Balay #undef __FUNCT__ 6304a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 631dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 632289bc588SBarry Smith { 633c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 63487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 635dfbe8321SBarry Smith PetscErrorCode ierr; 6360805154bSBarry Smith PetscBLASInt m, n, _One=1; 6373a40ed3dSBarry Smith 6383a40ed3dSBarry Smith PetscFunctionBegin; 639c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 640c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 641d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 642600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6431ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6441ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6458b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 6461ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6471ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 648dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6493a40ed3dSBarry Smith PetscFunctionReturn(0); 650289bc588SBarry Smith } 6516ee01492SSatish Balay 6524a2ae208SSatish Balay #undef __FUNCT__ 6534a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 654dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 655289bc588SBarry Smith { 656c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 65787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 658dfbe8321SBarry Smith PetscErrorCode ierr; 6590805154bSBarry Smith PetscBLASInt m, n, _One=1; 66087828ca2SBarry Smith PetscScalar _DOne=1.0; 6613a40ed3dSBarry Smith 6623a40ed3dSBarry Smith PetscFunctionBegin; 663c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 664c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 665d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 666600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6671ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6681ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6698b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 6701ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6711ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 672dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6733a40ed3dSBarry Smith PetscFunctionReturn(0); 674289bc588SBarry Smith } 675289bc588SBarry Smith 676289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6774a2ae208SSatish Balay #undef __FUNCT__ 6784a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 67913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 680289bc588SBarry Smith { 681c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 68287828ca2SBarry Smith PetscScalar *v; 6836849ba73SBarry Smith PetscErrorCode ierr; 68413f74950SBarry Smith PetscInt i; 68567e560aaSBarry Smith 6863a40ed3dSBarry Smith PetscFunctionBegin; 687d0f46423SBarry Smith *ncols = A->cmap->n; 688289bc588SBarry Smith if (cols) { 689785e854fSJed Brown ierr = PetscMalloc1((A->cmap->n+1),cols);CHKERRQ(ierr); 690d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 691289bc588SBarry Smith } 692289bc588SBarry Smith if (vals) { 693785e854fSJed Brown ierr = PetscMalloc1((A->cmap->n+1),vals);CHKERRQ(ierr); 694289bc588SBarry Smith v = mat->v + row; 695d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 696289bc588SBarry Smith } 6973a40ed3dSBarry Smith PetscFunctionReturn(0); 698289bc588SBarry Smith } 6996ee01492SSatish Balay 7004a2ae208SSatish Balay #undef __FUNCT__ 7014a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 70213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 703289bc588SBarry Smith { 704dfbe8321SBarry Smith PetscErrorCode ierr; 7056e111a19SKarl Rupp 706606d414cSSatish Balay PetscFunctionBegin; 707606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 708606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 7093a40ed3dSBarry Smith PetscFunctionReturn(0); 710289bc588SBarry Smith } 711289bc588SBarry Smith /* ----------------------------------------------------------------*/ 7124a2ae208SSatish Balay #undef __FUNCT__ 7134a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 71413f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 715289bc588SBarry Smith { 716c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 71713f74950SBarry Smith PetscInt i,j,idx=0; 718d6dfbf8fSBarry Smith 7193a40ed3dSBarry Smith PetscFunctionBegin; 72071fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 721289bc588SBarry Smith if (!mat->roworiented) { 722dbb450caSBarry Smith if (addv == INSERT_VALUES) { 723289bc588SBarry Smith for (j=0; j<n; j++) { 724cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7252515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 726e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 72758804f6dSBarry Smith #endif 728289bc588SBarry Smith for (i=0; i<m; i++) { 729cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7302515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 731e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 73258804f6dSBarry Smith #endif 733cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 734289bc588SBarry Smith } 735289bc588SBarry Smith } 7363a40ed3dSBarry Smith } else { 737289bc588SBarry Smith for (j=0; j<n; j++) { 738cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7392515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 740e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 74158804f6dSBarry Smith #endif 742289bc588SBarry Smith for (i=0; i<m; i++) { 743cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7442515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 745e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 74658804f6dSBarry Smith #endif 747cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 748289bc588SBarry Smith } 749289bc588SBarry Smith } 750289bc588SBarry Smith } 7513a40ed3dSBarry Smith } else { 752dbb450caSBarry Smith if (addv == INSERT_VALUES) { 753e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 754cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7552515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 756e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 75758804f6dSBarry Smith #endif 758e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 759cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7602515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 761e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 76258804f6dSBarry Smith #endif 763cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 764e8d4e0b9SBarry Smith } 765e8d4e0b9SBarry Smith } 7663a40ed3dSBarry Smith } else { 767289bc588SBarry Smith for (i=0; i<m; i++) { 768cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7692515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 770e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 77158804f6dSBarry Smith #endif 772289bc588SBarry Smith for (j=0; j<n; j++) { 773cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7742515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 775e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 77658804f6dSBarry Smith #endif 777cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 778289bc588SBarry Smith } 779289bc588SBarry Smith } 780289bc588SBarry Smith } 781e8d4e0b9SBarry Smith } 7823a40ed3dSBarry Smith PetscFunctionReturn(0); 783289bc588SBarry Smith } 784e8d4e0b9SBarry Smith 7854a2ae208SSatish Balay #undef __FUNCT__ 7864a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 78713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 788ae80bb75SLois Curfman McInnes { 789ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 79013f74950SBarry Smith PetscInt i,j; 791ae80bb75SLois Curfman McInnes 7923a40ed3dSBarry Smith PetscFunctionBegin; 793ae80bb75SLois Curfman McInnes /* row-oriented output */ 794ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 79597e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 796e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n); 797ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7986f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 799e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n); 80097e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 801ae80bb75SLois Curfman McInnes } 802ae80bb75SLois Curfman McInnes } 8033a40ed3dSBarry Smith PetscFunctionReturn(0); 804ae80bb75SLois Curfman McInnes } 805ae80bb75SLois Curfman McInnes 806289bc588SBarry Smith /* -----------------------------------------------------------------*/ 807289bc588SBarry Smith 8084a2ae208SSatish Balay #undef __FUNCT__ 8095bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 810112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 811aabbc4fbSShri Abhyankar { 812aabbc4fbSShri Abhyankar Mat_SeqDense *a; 813aabbc4fbSShri Abhyankar PetscErrorCode ierr; 814aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 815aabbc4fbSShri Abhyankar int fd; 816aabbc4fbSShri Abhyankar PetscMPIInt size; 817aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 818aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 819ce94432eSBarry Smith MPI_Comm comm; 820aabbc4fbSShri Abhyankar 821aabbc4fbSShri Abhyankar PetscFunctionBegin; 822ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 823aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 824aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 825aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 826aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 827aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 828aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 829aabbc4fbSShri Abhyankar 830aabbc4fbSShri Abhyankar /* set global size if not set already*/ 831aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 832aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 833aabbc4fbSShri Abhyankar } else { 834aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 835aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 836aabbc4fbSShri Abhyankar if (M != grows || N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols); 837aabbc4fbSShri Abhyankar } 8380298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 839aabbc4fbSShri Abhyankar 840aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 841aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 842aabbc4fbSShri Abhyankar v = a->v; 843aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 844aabbc4fbSShri Abhyankar from row major to column major */ 845785e854fSJed Brown ierr = PetscMalloc1((M*N > 0 ? M*N : 1),&w);CHKERRQ(ierr); 846aabbc4fbSShri Abhyankar /* read in nonzero values */ 847aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 848aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 849aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 850aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 851aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 852aabbc4fbSShri Abhyankar } 853aabbc4fbSShri Abhyankar } 854aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 855aabbc4fbSShri Abhyankar } else { 856aabbc4fbSShri Abhyankar /* read row lengths */ 857785e854fSJed Brown ierr = PetscMalloc1((M+1),&rowlengths);CHKERRQ(ierr); 858aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 859aabbc4fbSShri Abhyankar 860aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 861aabbc4fbSShri Abhyankar v = a->v; 862aabbc4fbSShri Abhyankar 863aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 864785e854fSJed Brown ierr = PetscMalloc1((nz+1),&scols);CHKERRQ(ierr); 865aabbc4fbSShri Abhyankar cols = scols; 866aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 867785e854fSJed Brown ierr = PetscMalloc1((nz+1),&svals);CHKERRQ(ierr); 868aabbc4fbSShri Abhyankar vals = svals; 869aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 870aabbc4fbSShri Abhyankar 871aabbc4fbSShri Abhyankar /* insert into matrix */ 872aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 873aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 874aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 875aabbc4fbSShri Abhyankar } 876aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 877aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 878aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 879aabbc4fbSShri Abhyankar } 880aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 881aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 882aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 883aabbc4fbSShri Abhyankar } 884aabbc4fbSShri Abhyankar 885aabbc4fbSShri Abhyankar #undef __FUNCT__ 8864a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8876849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 888289bc588SBarry Smith { 889932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 890dfbe8321SBarry Smith PetscErrorCode ierr; 89113f74950SBarry Smith PetscInt i,j; 8922dcb1b2aSMatthew Knepley const char *name; 89387828ca2SBarry Smith PetscScalar *v; 894f3ef73ceSBarry Smith PetscViewerFormat format; 8955f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 896ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8975f481a85SSatish Balay #endif 898932b0c3eSLois Curfman McInnes 8993a40ed3dSBarry Smith PetscFunctionBegin; 900b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 901456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 9023a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 903fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 904d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 905dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer);CHKERRQ(ierr); 906d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 90744cd7ae7SLois Curfman McInnes v = a->v + i; 90877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 909d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 910aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 911329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 91257622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr); 913329f5518SBarry Smith } else if (PetscRealPart(*v)) { 91457622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)PetscRealPart(*v));CHKERRQ(ierr); 9156831982aSBarry Smith } 91680cd9d93SLois Curfman McInnes #else 9176831982aSBarry Smith if (*v) { 91857622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",j,(double)*v);CHKERRQ(ierr); 9196831982aSBarry Smith } 92080cd9d93SLois Curfman McInnes #endif 9211b807ce4Svictorle v += a->lda; 92280cd9d93SLois Curfman McInnes } 923b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 92480cd9d93SLois Curfman McInnes } 925d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 9263a40ed3dSBarry Smith } else { 927d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 928aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 92947989497SBarry Smith /* determine if matrix has all real values */ 93047989497SBarry Smith v = a->v; 931d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 932ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;} 93347989497SBarry Smith } 93447989497SBarry Smith #endif 935fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 9363a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 937d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 938d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 939fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9407566de4bSShri Abhyankar } else { 941dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer);CHKERRQ(ierr); 942ffac6cdbSBarry Smith } 943ffac6cdbSBarry Smith 944d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 945932b0c3eSLois Curfman McInnes v = a->v + i; 946d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 947aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 94847989497SBarry Smith if (allreal) { 949c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));CHKERRQ(ierr); 95047989497SBarry Smith } else { 951c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr); 95247989497SBarry Smith } 953289bc588SBarry Smith #else 954c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);CHKERRQ(ierr); 955289bc588SBarry Smith #endif 9561b807ce4Svictorle v += a->lda; 957289bc588SBarry Smith } 958b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 959289bc588SBarry Smith } 960fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 961b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 962ffac6cdbSBarry Smith } 963d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 964da3a660dSBarry Smith } 965b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9663a40ed3dSBarry Smith PetscFunctionReturn(0); 967289bc588SBarry Smith } 968289bc588SBarry Smith 9694a2ae208SSatish Balay #undef __FUNCT__ 9704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9716849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 972932b0c3eSLois Curfman McInnes { 973932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9746849ba73SBarry Smith PetscErrorCode ierr; 97513f74950SBarry Smith int fd; 976d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 977f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 978f4403165SShri Abhyankar PetscViewerFormat format; 979932b0c3eSLois Curfman McInnes 9803a40ed3dSBarry Smith PetscFunctionBegin; 981b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 98290ace30eSBarry Smith 983f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 984f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 985f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 986785e854fSJed Brown ierr = PetscMalloc1(4,&col_lens);CHKERRQ(ierr); 9872205254eSKarl Rupp 988f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 989f4403165SShri Abhyankar col_lens[1] = m; 990f4403165SShri Abhyankar col_lens[2] = n; 991f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 9922205254eSKarl Rupp 993f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 994f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 995f4403165SShri Abhyankar 996f4403165SShri Abhyankar /* write out matrix, by rows */ 997785e854fSJed Brown ierr = PetscMalloc1((m*n+1),&vals);CHKERRQ(ierr); 998f4403165SShri Abhyankar v = a->v; 999f4403165SShri Abhyankar for (j=0; j<n; j++) { 1000f4403165SShri Abhyankar for (i=0; i<m; i++) { 1001f4403165SShri Abhyankar vals[j + i*n] = *v++; 1002f4403165SShri Abhyankar } 1003f4403165SShri Abhyankar } 1004f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1005f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 1006f4403165SShri Abhyankar } else { 1007785e854fSJed Brown ierr = PetscMalloc1((4+nz),&col_lens);CHKERRQ(ierr); 10082205254eSKarl Rupp 10090700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 1010932b0c3eSLois Curfman McInnes col_lens[1] = m; 1011932b0c3eSLois Curfman McInnes col_lens[2] = n; 1012932b0c3eSLois Curfman McInnes col_lens[3] = nz; 1013932b0c3eSLois Curfman McInnes 1014932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 1015932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 10166f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1017932b0c3eSLois Curfman McInnes 1018932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 1019932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 1020932b0c3eSLois Curfman McInnes ict = 0; 1021932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1022932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 1023932b0c3eSLois Curfman McInnes } 10246f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1025606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 1026932b0c3eSLois Curfman McInnes 1027932b0c3eSLois Curfman McInnes /* store nonzero values */ 1028785e854fSJed Brown ierr = PetscMalloc1((nz+1),&anonz);CHKERRQ(ierr); 1029932b0c3eSLois Curfman McInnes ict = 0; 1030932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1031932b0c3eSLois Curfman McInnes v = a->v + i; 1032932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 10331b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 1034932b0c3eSLois Curfman McInnes } 1035932b0c3eSLois Curfman McInnes } 10366f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1037606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 1038f4403165SShri Abhyankar } 10393a40ed3dSBarry Smith PetscFunctionReturn(0); 1040932b0c3eSLois Curfman McInnes } 1041932b0c3eSLois Curfman McInnes 10429804daf3SBarry Smith #include <petscdraw.h> 10434a2ae208SSatish Balay #undef __FUNCT__ 10444a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1045dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1046f1af5d2fSBarry Smith { 1047f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1048f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10496849ba73SBarry Smith PetscErrorCode ierr; 1050d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 105187828ca2SBarry Smith PetscScalar *v = a->v; 1052b0a32e0cSBarry Smith PetscViewer viewer; 1053b0a32e0cSBarry Smith PetscDraw popup; 1054329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1055f3ef73ceSBarry Smith PetscViewerFormat format; 1056f1af5d2fSBarry Smith 1057f1af5d2fSBarry Smith PetscFunctionBegin; 1058f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1059b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1060b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1061f1af5d2fSBarry Smith 1062f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1063fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1064f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1065b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1066f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1067f1af5d2fSBarry Smith x_l = j; 1068f1af5d2fSBarry Smith x_r = x_l + 1.0; 1069f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1070f1af5d2fSBarry Smith y_l = m - i - 1.0; 1071f1af5d2fSBarry Smith y_r = y_l + 1.0; 1072329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1073b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1074329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1075b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1076f1af5d2fSBarry Smith } else { 1077f1af5d2fSBarry Smith continue; 1078f1af5d2fSBarry Smith } 1079b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1080f1af5d2fSBarry Smith } 1081f1af5d2fSBarry Smith } 1082f1af5d2fSBarry Smith } else { 1083f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1084f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1085f1af5d2fSBarry Smith for (i = 0; i < m*n; i++) { 1086f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1087f1af5d2fSBarry Smith } 1088b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1089b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1090b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1091f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1092f1af5d2fSBarry Smith x_l = j; 1093f1af5d2fSBarry Smith x_r = x_l + 1.0; 1094f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1095f1af5d2fSBarry Smith y_l = m - i - 1.0; 1096f1af5d2fSBarry Smith y_r = y_l + 1.0; 1097b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1098b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1099f1af5d2fSBarry Smith } 1100f1af5d2fSBarry Smith } 1101f1af5d2fSBarry Smith } 1102f1af5d2fSBarry Smith PetscFunctionReturn(0); 1103f1af5d2fSBarry Smith } 1104f1af5d2fSBarry Smith 11054a2ae208SSatish Balay #undef __FUNCT__ 11064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1107dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1108f1af5d2fSBarry Smith { 1109b0a32e0cSBarry Smith PetscDraw draw; 1110ace3abfcSBarry Smith PetscBool isnull; 1111329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1112dfbe8321SBarry Smith PetscErrorCode ierr; 1113f1af5d2fSBarry Smith 1114f1af5d2fSBarry Smith PetscFunctionBegin; 1115b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1116b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1117abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1118f1af5d2fSBarry Smith 1119f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1120d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1121f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1122b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1123b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 11240298fd71SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr); 1125f1af5d2fSBarry Smith PetscFunctionReturn(0); 1126f1af5d2fSBarry Smith } 1127f1af5d2fSBarry Smith 11284a2ae208SSatish Balay #undef __FUNCT__ 11294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1130dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1131932b0c3eSLois Curfman McInnes { 1132dfbe8321SBarry Smith PetscErrorCode ierr; 1133ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1134932b0c3eSLois Curfman McInnes 11353a40ed3dSBarry Smith PetscFunctionBegin; 1136251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1137251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1138251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11390f5bd95cSBarry Smith 1140c45a1595SBarry Smith if (iascii) { 1141c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11420f5bd95cSBarry Smith } else if (isbinary) { 11433a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1144f1af5d2fSBarry Smith } else if (isdraw) { 1145f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 1146932b0c3eSLois Curfman McInnes } 11473a40ed3dSBarry Smith PetscFunctionReturn(0); 1148932b0c3eSLois Curfman McInnes } 1149289bc588SBarry Smith 11504a2ae208SSatish Balay #undef __FUNCT__ 11514a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1152dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1153289bc588SBarry Smith { 1154ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1155dfbe8321SBarry Smith PetscErrorCode ierr; 115690f02eecSBarry Smith 11573a40ed3dSBarry Smith PetscFunctionBegin; 1158aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1159d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1160a5a9c739SBarry Smith #endif 116105b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11626857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1163bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1164dbd8c25aSHong Zhang 1165dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1166bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL);CHKERRQ(ierr); 1167bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL);CHKERRQ(ierr); 1168bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL);CHKERRQ(ierr); 1169bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 1170bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 1171bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 11723bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMult_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 11733bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 11743bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)mat,"MatTransposeMatMultNumeric_seqaij_seqdense_C",NULL);CHKERRQ(ierr); 11753a40ed3dSBarry Smith PetscFunctionReturn(0); 1176289bc588SBarry Smith } 1177289bc588SBarry Smith 11784a2ae208SSatish Balay #undef __FUNCT__ 11794a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1180fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1181289bc588SBarry Smith { 1182c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11836849ba73SBarry Smith PetscErrorCode ierr; 118413f74950SBarry Smith PetscInt k,j,m,n,M; 118587828ca2SBarry Smith PetscScalar *v,tmp; 118648b35521SBarry Smith 11873a40ed3dSBarry Smith PetscFunctionBegin; 1188d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1189e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1190e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1191e7e72b3dSBarry Smith else { 1192d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1193289bc588SBarry Smith for (k=0; k<j; k++) { 11941b807ce4Svictorle tmp = v[j + k*M]; 11951b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11961b807ce4Svictorle v[k + j*M] = tmp; 1197289bc588SBarry Smith } 1198289bc588SBarry Smith } 1199d64ed03dSBarry Smith } 12003a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1201d3e5ee88SLois Curfman McInnes Mat tmat; 1202ec8511deSBarry Smith Mat_SeqDense *tmatd; 120387828ca2SBarry Smith PetscScalar *v2; 1204ea709b57SSatish Balay 1205fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 1206ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&tmat);CHKERRQ(ierr); 1207d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 12087adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 12090298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(tmat,NULL);CHKERRQ(ierr); 1210fc4dec0aSBarry Smith } else { 1211fc4dec0aSBarry Smith tmat = *matout; 1212fc4dec0aSBarry Smith } 1213ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 12140de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1215d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 12161b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1217d3e5ee88SLois Curfman McInnes } 12186d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12196d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1220d3e5ee88SLois Curfman McInnes *matout = tmat; 122148b35521SBarry Smith } 12223a40ed3dSBarry Smith PetscFunctionReturn(0); 1223289bc588SBarry Smith } 1224289bc588SBarry Smith 12254a2ae208SSatish Balay #undef __FUNCT__ 12264a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1227ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1228289bc588SBarry Smith { 1229c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1230c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 123113f74950SBarry Smith PetscInt i,j; 1232a2ea699eSBarry Smith PetscScalar *v1,*v2; 12339ea5d5aeSSatish Balay 12343a40ed3dSBarry Smith PetscFunctionBegin; 1235d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1236d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1237d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12381b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1239d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12403a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12411b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12421b807ce4Svictorle } 1243289bc588SBarry Smith } 124477c4ece6SBarry Smith *flg = PETSC_TRUE; 12453a40ed3dSBarry Smith PetscFunctionReturn(0); 1246289bc588SBarry Smith } 1247289bc588SBarry Smith 12484a2ae208SSatish Balay #undef __FUNCT__ 12494a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1250dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1251289bc588SBarry Smith { 1252c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1253dfbe8321SBarry Smith PetscErrorCode ierr; 125413f74950SBarry Smith PetscInt i,n,len; 125587828ca2SBarry Smith PetscScalar *x,zero = 0.0; 125644cd7ae7SLois Curfman McInnes 12573a40ed3dSBarry Smith PetscFunctionBegin; 12582dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12597a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12601ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1261d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1262e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 126344cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12641b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1265289bc588SBarry Smith } 12661ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12673a40ed3dSBarry Smith PetscFunctionReturn(0); 1268289bc588SBarry Smith } 1269289bc588SBarry Smith 12704a2ae208SSatish Balay #undef __FUNCT__ 12714a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1272dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1273289bc588SBarry Smith { 1274c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 127587828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1276dfbe8321SBarry Smith PetscErrorCode ierr; 1277d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 127855659b69SBarry Smith 12793a40ed3dSBarry Smith PetscFunctionBegin; 128028988994SBarry Smith if (ll) { 12817a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12821ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1283e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1284da3a660dSBarry Smith for (i=0; i<m; i++) { 1285da3a660dSBarry Smith x = l[i]; 1286da3a660dSBarry Smith v = mat->v + i; 1287da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1288da3a660dSBarry Smith } 12891ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1290efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1291da3a660dSBarry Smith } 129228988994SBarry Smith if (rr) { 12937a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12941ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1295e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1296da3a660dSBarry Smith for (i=0; i<n; i++) { 1297da3a660dSBarry Smith x = r[i]; 1298da3a660dSBarry Smith v = mat->v + i*m; 12992205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 1300da3a660dSBarry Smith } 13011ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1302efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1303da3a660dSBarry Smith } 13043a40ed3dSBarry Smith PetscFunctionReturn(0); 1305289bc588SBarry Smith } 1306289bc588SBarry Smith 13074a2ae208SSatish Balay #undef __FUNCT__ 13084a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1309dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1310289bc588SBarry Smith { 1311c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 131287828ca2SBarry Smith PetscScalar *v = mat->v; 1313329f5518SBarry Smith PetscReal sum = 0.0; 1314d0f46423SBarry Smith PetscInt lda =mat->lda,m=A->rmap->n,i,j; 1315efee365bSSatish Balay PetscErrorCode ierr; 131655659b69SBarry Smith 13173a40ed3dSBarry Smith PetscFunctionBegin; 1318289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1319a5ce6ee0Svictorle if (lda>m) { 1320d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1321a5ce6ee0Svictorle v = mat->v+j*lda; 1322a5ce6ee0Svictorle for (i=0; i<m; i++) { 1323a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1324a5ce6ee0Svictorle } 1325a5ce6ee0Svictorle } 1326a5ce6ee0Svictorle } else { 1327d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1328329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1329289bc588SBarry Smith } 1330a5ce6ee0Svictorle } 13318f1a2a5eSBarry Smith *nrm = PetscSqrtReal(sum); 1332dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13333a40ed3dSBarry Smith } else if (type == NORM_1) { 1334064f8208SBarry Smith *nrm = 0.0; 1335d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13361b807ce4Svictorle v = mat->v + j*mat->lda; 1337289bc588SBarry Smith sum = 0.0; 1338d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 133933a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1340289bc588SBarry Smith } 1341064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1342289bc588SBarry Smith } 1343d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13443a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1345064f8208SBarry Smith *nrm = 0.0; 1346d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1347289bc588SBarry Smith v = mat->v + j; 1348289bc588SBarry Smith sum = 0.0; 1349d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13501b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1351289bc588SBarry Smith } 1352064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1353289bc588SBarry Smith } 1354d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1355e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13563a40ed3dSBarry Smith PetscFunctionReturn(0); 1357289bc588SBarry Smith } 1358289bc588SBarry Smith 13594a2ae208SSatish Balay #undef __FUNCT__ 13604a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1361ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1362289bc588SBarry Smith { 1363c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 136463ba0a88SBarry Smith PetscErrorCode ierr; 136567e560aaSBarry Smith 13663a40ed3dSBarry Smith PetscFunctionBegin; 1367b5a2b587SKris Buschelman switch (op) { 1368b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13694e0d8c25SBarry Smith aij->roworiented = flg; 1370b5a2b587SKris Buschelman break; 1371512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1372b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13733971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13744e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 137513fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 1376b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1377b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 13785021d80fSJed Brown case MAT_IGNORE_LOWER_TRIANGULAR: 13795021d80fSJed Brown ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 13805021d80fSJed Brown break; 13815021d80fSJed Brown case MAT_SPD: 138277e54ba9SKris Buschelman case MAT_SYMMETRIC: 138377e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13849a4540c5SBarry Smith case MAT_HERMITIAN: 13859a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 13865021d80fSJed Brown /* These options are handled directly by MatSetOption() */ 138777e54ba9SKris Buschelman break; 1388b5a2b587SKris Buschelman default: 1389e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13903a40ed3dSBarry Smith } 13913a40ed3dSBarry Smith PetscFunctionReturn(0); 1392289bc588SBarry Smith } 1393289bc588SBarry Smith 13944a2ae208SSatish Balay #undef __FUNCT__ 13954a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1396dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13976f0a148fSBarry Smith { 1398ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13996849ba73SBarry Smith PetscErrorCode ierr; 1400d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 14013a40ed3dSBarry Smith 14023a40ed3dSBarry Smith PetscFunctionBegin; 1403a5ce6ee0Svictorle if (lda>m) { 1404d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1405a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1406a5ce6ee0Svictorle } 1407a5ce6ee0Svictorle } else { 1408d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1409a5ce6ee0Svictorle } 14103a40ed3dSBarry Smith PetscFunctionReturn(0); 14116f0a148fSBarry Smith } 14126f0a148fSBarry Smith 14134a2ae208SSatish Balay #undef __FUNCT__ 14144a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 14152b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 14166f0a148fSBarry Smith { 141797b48c8fSBarry Smith PetscErrorCode ierr; 1418ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1419b9679d65SBarry Smith PetscInt m = l->lda, n = A->cmap->n, i,j; 142097b48c8fSBarry Smith PetscScalar *slot,*bb; 142197b48c8fSBarry Smith const PetscScalar *xx; 142255659b69SBarry Smith 14233a40ed3dSBarry Smith PetscFunctionBegin; 1424b9679d65SBarry Smith #if defined(PETSC_USE_DEBUG) 1425b9679d65SBarry Smith for (i=0; i<N; i++) { 1426b9679d65SBarry Smith if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 1427b9679d65SBarry Smith if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n); 1428b9679d65SBarry Smith } 1429b9679d65SBarry Smith #endif 1430b9679d65SBarry Smith 143197b48c8fSBarry Smith /* fix right hand side if needed */ 143297b48c8fSBarry Smith if (x && b) { 143397b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 143497b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 14352205254eSKarl Rupp for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]]; 143697b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 143797b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 143897b48c8fSBarry Smith } 143997b48c8fSBarry Smith 14406f0a148fSBarry Smith for (i=0; i<N; i++) { 14416f0a148fSBarry Smith slot = l->v + rows[i]; 1442b9679d65SBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 14436f0a148fSBarry Smith } 1444f4df32b1SMatthew Knepley if (diag != 0.0) { 1445b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 14466f0a148fSBarry Smith for (i=0; i<N; i++) { 1447b9679d65SBarry Smith slot = l->v + (m+1)*rows[i]; 1448f4df32b1SMatthew Knepley *slot = diag; 14496f0a148fSBarry Smith } 14506f0a148fSBarry Smith } 14513a40ed3dSBarry Smith PetscFunctionReturn(0); 14526f0a148fSBarry Smith } 1453557bce09SLois Curfman McInnes 14544a2ae208SSatish Balay #undef __FUNCT__ 14558c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_SeqDense" 14568c778c55SBarry Smith PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[]) 145764e87e97SBarry Smith { 1458c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14593a40ed3dSBarry Smith 14603a40ed3dSBarry Smith PetscFunctionBegin; 1461e32f2f54SBarry Smith if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 146264e87e97SBarry Smith *array = mat->v; 14633a40ed3dSBarry Smith PetscFunctionReturn(0); 146464e87e97SBarry Smith } 14650754003eSLois Curfman McInnes 14664a2ae208SSatish Balay #undef __FUNCT__ 14678c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_SeqDense" 14688c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1469ff14e315SSatish Balay { 14703a40ed3dSBarry Smith PetscFunctionBegin; 147109b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14723a40ed3dSBarry Smith PetscFunctionReturn(0); 1473ff14e315SSatish Balay } 14740754003eSLois Curfman McInnes 14754a2ae208SSatish Balay #undef __FUNCT__ 14768c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray" 1477dec5eb66SMatthew G Knepley /*@C 14788c778c55SBarry Smith MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored 147973a71a0fSBarry Smith 148073a71a0fSBarry Smith Not Collective 148173a71a0fSBarry Smith 148273a71a0fSBarry Smith Input Parameter: 148373a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 148473a71a0fSBarry Smith 148573a71a0fSBarry Smith Output Parameter: 148673a71a0fSBarry Smith . array - pointer to the data 148773a71a0fSBarry Smith 148873a71a0fSBarry Smith Level: intermediate 148973a71a0fSBarry Smith 14908c778c55SBarry Smith .seealso: MatDenseRestoreArray() 149173a71a0fSBarry Smith @*/ 14928c778c55SBarry Smith PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array) 149373a71a0fSBarry Smith { 149473a71a0fSBarry Smith PetscErrorCode ierr; 149573a71a0fSBarry Smith 149673a71a0fSBarry Smith PetscFunctionBegin; 14978c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 149873a71a0fSBarry Smith PetscFunctionReturn(0); 149973a71a0fSBarry Smith } 150073a71a0fSBarry Smith 150173a71a0fSBarry Smith #undef __FUNCT__ 15028c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray" 1503dec5eb66SMatthew G Knepley /*@C 15048c778c55SBarry Smith MatDenseRestoreArray - returns access to the array where the data for a SeqDense matrix is stored obtained by MatDenseGetArray() 150573a71a0fSBarry Smith 150673a71a0fSBarry Smith Not Collective 150773a71a0fSBarry Smith 150873a71a0fSBarry Smith Input Parameters: 150973a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 151073a71a0fSBarry Smith . array - pointer to the data 151173a71a0fSBarry Smith 151273a71a0fSBarry Smith Level: intermediate 151373a71a0fSBarry Smith 15148c778c55SBarry Smith .seealso: MatDenseGetArray() 151573a71a0fSBarry Smith @*/ 15168c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array) 151773a71a0fSBarry Smith { 151873a71a0fSBarry Smith PetscErrorCode ierr; 151973a71a0fSBarry Smith 152073a71a0fSBarry Smith PetscFunctionBegin; 15218c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 152273a71a0fSBarry Smith PetscFunctionReturn(0); 152373a71a0fSBarry Smith } 152473a71a0fSBarry Smith 152573a71a0fSBarry Smith #undef __FUNCT__ 15264a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 152713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 15280754003eSLois Curfman McInnes { 1529c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 15306849ba73SBarry Smith PetscErrorCode ierr; 15315d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 15325d0c19d7SBarry Smith const PetscInt *irow,*icol; 153387828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 15340754003eSLois Curfman McInnes Mat newmat; 15350754003eSLois Curfman McInnes 15363a40ed3dSBarry Smith PetscFunctionBegin; 153778b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 153878b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1539e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1540e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 15410754003eSLois Curfman McInnes 1542182d2002SSatish Balay /* Check submatrixcall */ 1543182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 154413f74950SBarry Smith PetscInt n_cols,n_rows; 1545182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 154621a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1547f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1548c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 154921a2c019SBarry Smith } 1550182d2002SSatish Balay newmat = *B; 1551182d2002SSatish Balay } else { 15520754003eSLois Curfman McInnes /* Create and fill new matrix */ 1553ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr); 1554f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 15557adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 15560298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 1557182d2002SSatish Balay } 1558182d2002SSatish Balay 1559182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1560182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1561182d2002SSatish Balay 1562182d2002SSatish Balay for (i=0; i<ncols; i++) { 15636de62eeeSBarry Smith av = v + mat->lda*icol[i]; 15642205254eSKarl Rupp for (j=0; j<nrows; j++) *bv++ = av[irow[j]]; 15650754003eSLois Curfman McInnes } 1566182d2002SSatish Balay 1567182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 15686d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15696d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15700754003eSLois Curfman McInnes 15710754003eSLois Curfman McInnes /* Free work space */ 157278b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 157378b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1574182d2002SSatish Balay *B = newmat; 15753a40ed3dSBarry Smith PetscFunctionReturn(0); 15760754003eSLois Curfman McInnes } 15770754003eSLois Curfman McInnes 15784a2ae208SSatish Balay #undef __FUNCT__ 15794a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 158013f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1581905e6a2fSBarry Smith { 15826849ba73SBarry Smith PetscErrorCode ierr; 158313f74950SBarry Smith PetscInt i; 1584905e6a2fSBarry Smith 15853a40ed3dSBarry Smith PetscFunctionBegin; 1586905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1587785e854fSJed Brown ierr = PetscMalloc1((n+1),B);CHKERRQ(ierr); 1588905e6a2fSBarry Smith } 1589905e6a2fSBarry Smith 1590905e6a2fSBarry Smith for (i=0; i<n; i++) { 15916a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1592905e6a2fSBarry Smith } 15933a40ed3dSBarry Smith PetscFunctionReturn(0); 1594905e6a2fSBarry Smith } 1595905e6a2fSBarry Smith 15964a2ae208SSatish Balay #undef __FUNCT__ 1597c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1598c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1599c0aa2d19SHong Zhang { 1600c0aa2d19SHong Zhang PetscFunctionBegin; 1601c0aa2d19SHong Zhang PetscFunctionReturn(0); 1602c0aa2d19SHong Zhang } 1603c0aa2d19SHong Zhang 1604c0aa2d19SHong Zhang #undef __FUNCT__ 1605c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1606c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1607c0aa2d19SHong Zhang { 1608c0aa2d19SHong Zhang PetscFunctionBegin; 1609c0aa2d19SHong Zhang PetscFunctionReturn(0); 1610c0aa2d19SHong Zhang } 1611c0aa2d19SHong Zhang 1612c0aa2d19SHong Zhang #undef __FUNCT__ 16134a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1614dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 16154b0e389bSBarry Smith { 16164b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data; 16176849ba73SBarry Smith PetscErrorCode ierr; 1618d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 16193a40ed3dSBarry Smith 16203a40ed3dSBarry Smith PetscFunctionBegin; 162133f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 162233f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1623cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 16243a40ed3dSBarry Smith PetscFunctionReturn(0); 16253a40ed3dSBarry Smith } 1626e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1627a5ce6ee0Svictorle if (lda1>m || lda2>m) { 16280dbb7854Svictorle for (j=0; j<n; j++) { 1629a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1630a5ce6ee0Svictorle } 1631a5ce6ee0Svictorle } else { 1632d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1633a5ce6ee0Svictorle } 1634273d9f13SBarry Smith PetscFunctionReturn(0); 1635273d9f13SBarry Smith } 1636273d9f13SBarry Smith 16374a2ae208SSatish Balay #undef __FUNCT__ 16384994cf47SJed Brown #define __FUNCT__ "MatSetUp_SeqDense" 16394994cf47SJed Brown PetscErrorCode MatSetUp_SeqDense(Mat A) 1640273d9f13SBarry Smith { 1641dfbe8321SBarry Smith PetscErrorCode ierr; 1642273d9f13SBarry Smith 1643273d9f13SBarry Smith PetscFunctionBegin; 1644273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 16453a40ed3dSBarry Smith PetscFunctionReturn(0); 16464b0e389bSBarry Smith } 16474b0e389bSBarry Smith 1648284134d9SBarry Smith #undef __FUNCT__ 1649ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1650ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1651ba337c44SJed Brown { 1652ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1653ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1654ba337c44SJed Brown PetscScalar *aa = a->v; 1655ba337c44SJed Brown 1656ba337c44SJed Brown PetscFunctionBegin; 1657ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1658ba337c44SJed Brown PetscFunctionReturn(0); 1659ba337c44SJed Brown } 1660ba337c44SJed Brown 1661ba337c44SJed Brown #undef __FUNCT__ 1662ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1663ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1664ba337c44SJed Brown { 1665ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1666ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1667ba337c44SJed Brown PetscScalar *aa = a->v; 1668ba337c44SJed Brown 1669ba337c44SJed Brown PetscFunctionBegin; 1670ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1671ba337c44SJed Brown PetscFunctionReturn(0); 1672ba337c44SJed Brown } 1673ba337c44SJed Brown 1674ba337c44SJed Brown #undef __FUNCT__ 1675ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1676ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1677ba337c44SJed Brown { 1678ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1679ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1680ba337c44SJed Brown PetscScalar *aa = a->v; 1681ba337c44SJed Brown 1682ba337c44SJed Brown PetscFunctionBegin; 1683ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1684ba337c44SJed Brown PetscFunctionReturn(0); 1685ba337c44SJed Brown } 1686284134d9SBarry Smith 1687a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1688a9fe9ddaSSatish Balay #undef __FUNCT__ 1689a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1690a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1691a9fe9ddaSSatish Balay { 1692a9fe9ddaSSatish Balay PetscErrorCode ierr; 1693a9fe9ddaSSatish Balay 1694a9fe9ddaSSatish Balay PetscFunctionBegin; 1695a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX) { 16963ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1697a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 16983ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1699a9fe9ddaSSatish Balay } 17003ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1701a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 17023ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1703a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1704a9fe9ddaSSatish Balay } 1705a9fe9ddaSSatish Balay 1706a9fe9ddaSSatish Balay #undef __FUNCT__ 1707a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1708a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1709a9fe9ddaSSatish Balay { 1710ee16a9a1SHong Zhang PetscErrorCode ierr; 1711d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1712ee16a9a1SHong Zhang Mat Cmat; 1713a9fe9ddaSSatish Balay 1714ee16a9a1SHong Zhang PetscFunctionBegin; 1715e32f2f54SBarry Smith if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n); 1716ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1717ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1718ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 17190298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr); 1720d73949e8SHong Zhang 1721ee16a9a1SHong Zhang *C = Cmat; 1722ee16a9a1SHong Zhang PetscFunctionReturn(0); 1723ee16a9a1SHong Zhang } 1724a9fe9ddaSSatish Balay 172598a3b096SSatish Balay #undef __FUNCT__ 1726a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1727a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1728a9fe9ddaSSatish Balay { 1729a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1730a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1731a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17320805154bSBarry Smith PetscBLASInt m,n,k; 1733a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1734c5df96a5SBarry Smith PetscErrorCode ierr; 1735a9fe9ddaSSatish Balay 1736a9fe9ddaSSatish Balay PetscFunctionBegin; 1737c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 1738c5df96a5SBarry Smith ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr); 1739c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&k);CHKERRQ(ierr); 17408b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda)); 1741a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1742a9fe9ddaSSatish Balay } 1743a9fe9ddaSSatish Balay 1744a9fe9ddaSSatish Balay #undef __FUNCT__ 174575648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqDense_SeqDense" 174675648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1747a9fe9ddaSSatish Balay { 1748a9fe9ddaSSatish Balay PetscErrorCode ierr; 1749a9fe9ddaSSatish Balay 1750a9fe9ddaSSatish Balay PetscFunctionBegin; 1751a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX) { 17523ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 175375648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 17543ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1755a9fe9ddaSSatish Balay } 17563ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 175775648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 17583ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1759a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1760a9fe9ddaSSatish Balay } 1761a9fe9ddaSSatish Balay 1762a9fe9ddaSSatish Balay #undef __FUNCT__ 176375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqDense_SeqDense" 176475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1765a9fe9ddaSSatish Balay { 1766ee16a9a1SHong Zhang PetscErrorCode ierr; 1767d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1768ee16a9a1SHong Zhang Mat Cmat; 1769a9fe9ddaSSatish Balay 1770ee16a9a1SHong Zhang PetscFunctionBegin; 1771e32f2f54SBarry Smith if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n); 1772ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1773ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1774ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 17750298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr); 17762205254eSKarl Rupp 1777ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 17782205254eSKarl Rupp 1779ee16a9a1SHong Zhang *C = Cmat; 1780ee16a9a1SHong Zhang PetscFunctionReturn(0); 1781ee16a9a1SHong Zhang } 1782a9fe9ddaSSatish Balay 1783a9fe9ddaSSatish Balay #undef __FUNCT__ 178475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqDense_SeqDense" 178575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1786a9fe9ddaSSatish Balay { 1787a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1788a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1789a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17900805154bSBarry Smith PetscBLASInt m,n,k; 1791a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1792c5df96a5SBarry Smith PetscErrorCode ierr; 1793a9fe9ddaSSatish Balay 1794a9fe9ddaSSatish Balay PetscFunctionBegin; 1795c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&m);CHKERRQ(ierr); 1796c5df96a5SBarry Smith ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr); 1797c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&k);CHKERRQ(ierr); 17982fbe02b9SBarry Smith /* 17992fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 18002fbe02b9SBarry Smith */ 18018b83055fSJed Brown PetscStackCallBLAS("BLASgemv",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda)); 1802a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1803a9fe9ddaSSatish Balay } 1804985db425SBarry Smith 1805985db425SBarry Smith #undef __FUNCT__ 1806985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1807985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1808985db425SBarry Smith { 1809985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1810985db425SBarry Smith PetscErrorCode ierr; 1811d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1812985db425SBarry Smith PetscScalar *x; 1813985db425SBarry Smith MatScalar *aa = a->v; 1814985db425SBarry Smith 1815985db425SBarry Smith PetscFunctionBegin; 1816e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1817985db425SBarry Smith 1818985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1819985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1820985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1821e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1822985db425SBarry Smith for (i=0; i<m; i++) { 1823985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1824985db425SBarry Smith for (j=1; j<n; j++) { 1825985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1826985db425SBarry Smith } 1827985db425SBarry Smith } 1828985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1829985db425SBarry Smith PetscFunctionReturn(0); 1830985db425SBarry Smith } 1831985db425SBarry Smith 1832985db425SBarry Smith #undef __FUNCT__ 1833985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1834985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1835985db425SBarry Smith { 1836985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1837985db425SBarry Smith PetscErrorCode ierr; 1838d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1839985db425SBarry Smith PetscScalar *x; 1840985db425SBarry Smith PetscReal atmp; 1841985db425SBarry Smith MatScalar *aa = a->v; 1842985db425SBarry Smith 1843985db425SBarry Smith PetscFunctionBegin; 1844e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1845985db425SBarry Smith 1846985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1847985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1848985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1849e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1850985db425SBarry Smith for (i=0; i<m; i++) { 18519189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1852985db425SBarry Smith for (j=1; j<n; j++) { 1853985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1854985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1855985db425SBarry Smith } 1856985db425SBarry Smith } 1857985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1858985db425SBarry Smith PetscFunctionReturn(0); 1859985db425SBarry Smith } 1860985db425SBarry Smith 1861985db425SBarry Smith #undef __FUNCT__ 1862985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1863985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1864985db425SBarry Smith { 1865985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1866985db425SBarry Smith PetscErrorCode ierr; 1867d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1868985db425SBarry Smith PetscScalar *x; 1869985db425SBarry Smith MatScalar *aa = a->v; 1870985db425SBarry Smith 1871985db425SBarry Smith PetscFunctionBegin; 1872e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1873985db425SBarry Smith 1874985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1875985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1876985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1877e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1878985db425SBarry Smith for (i=0; i<m; i++) { 1879985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1880985db425SBarry Smith for (j=1; j<n; j++) { 1881985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1882985db425SBarry Smith } 1883985db425SBarry Smith } 1884985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1885985db425SBarry Smith PetscFunctionReturn(0); 1886985db425SBarry Smith } 1887985db425SBarry Smith 18888d0534beSBarry Smith #undef __FUNCT__ 18898d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18908d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18918d0534beSBarry Smith { 18928d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18938d0534beSBarry Smith PetscErrorCode ierr; 18948d0534beSBarry Smith PetscScalar *x; 18958d0534beSBarry Smith 18968d0534beSBarry Smith PetscFunctionBegin; 1897e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18988d0534beSBarry Smith 18998d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1900d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 19018d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 19028d0534beSBarry Smith PetscFunctionReturn(0); 19038d0534beSBarry Smith } 19048d0534beSBarry Smith 19050716a85fSBarry Smith 19060716a85fSBarry Smith #undef __FUNCT__ 19070716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 19080716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 19090716a85fSBarry Smith { 19100716a85fSBarry Smith PetscErrorCode ierr; 19110716a85fSBarry Smith PetscInt i,j,m,n; 19120716a85fSBarry Smith PetscScalar *a; 19130716a85fSBarry Smith 19140716a85fSBarry Smith PetscFunctionBegin; 19150716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 19160716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 19178c778c55SBarry Smith ierr = MatDenseGetArray(A,&a);CHKERRQ(ierr); 19180716a85fSBarry Smith if (type == NORM_2) { 19190716a85fSBarry Smith for (i=0; i<n; i++) { 19200716a85fSBarry Smith for (j=0; j<m; j++) { 19210716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 19220716a85fSBarry Smith } 19230716a85fSBarry Smith a += m; 19240716a85fSBarry Smith } 19250716a85fSBarry Smith } else if (type == NORM_1) { 19260716a85fSBarry Smith for (i=0; i<n; i++) { 19270716a85fSBarry Smith for (j=0; j<m; j++) { 19280716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 19290716a85fSBarry Smith } 19300716a85fSBarry Smith a += m; 19310716a85fSBarry Smith } 19320716a85fSBarry Smith } else if (type == NORM_INFINITY) { 19330716a85fSBarry Smith for (i=0; i<n; i++) { 19340716a85fSBarry Smith for (j=0; j<m; j++) { 19350716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 19360716a85fSBarry Smith } 19370716a85fSBarry Smith a += m; 19380716a85fSBarry Smith } 1939ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType"); 19408c778c55SBarry Smith ierr = MatDenseRestoreArray(A,&a);CHKERRQ(ierr); 19410716a85fSBarry Smith if (type == NORM_2) { 19428f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 19430716a85fSBarry Smith } 19440716a85fSBarry Smith PetscFunctionReturn(0); 19450716a85fSBarry Smith } 19460716a85fSBarry Smith 194773a71a0fSBarry Smith #undef __FUNCT__ 194873a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_SeqDense" 194973a71a0fSBarry Smith static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx) 195073a71a0fSBarry Smith { 195173a71a0fSBarry Smith PetscErrorCode ierr; 195273a71a0fSBarry Smith PetscScalar *a; 195373a71a0fSBarry Smith PetscInt m,n,i; 195473a71a0fSBarry Smith 195573a71a0fSBarry Smith PetscFunctionBegin; 195673a71a0fSBarry Smith ierr = MatGetSize(x,&m,&n);CHKERRQ(ierr); 19578c778c55SBarry Smith ierr = MatDenseGetArray(x,&a);CHKERRQ(ierr); 195873a71a0fSBarry Smith for (i=0; i<m*n; i++) { 195973a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 196073a71a0fSBarry Smith } 19618c778c55SBarry Smith ierr = MatDenseRestoreArray(x,&a);CHKERRQ(ierr); 196273a71a0fSBarry Smith PetscFunctionReturn(0); 196373a71a0fSBarry Smith } 196473a71a0fSBarry Smith 196573a71a0fSBarry Smith 1966289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1967a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_SeqDense, 1968905e6a2fSBarry Smith MatGetRow_SeqDense, 1969905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1970905e6a2fSBarry Smith MatMult_SeqDense, 197197304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 19727c922b88SBarry Smith MatMultTranspose_SeqDense, 19737c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1974db4efbfdSBarry Smith 0, 1975db4efbfdSBarry Smith 0, 1976db4efbfdSBarry Smith 0, 1977db4efbfdSBarry Smith /* 10*/ 0, 1978905e6a2fSBarry Smith MatLUFactor_SeqDense, 1979905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 198041f059aeSBarry Smith MatSOR_SeqDense, 1981ec8511deSBarry Smith MatTranspose_SeqDense, 198297304618SKris Buschelman /* 15*/ MatGetInfo_SeqDense, 1983905e6a2fSBarry Smith MatEqual_SeqDense, 1984905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1985905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1986905e6a2fSBarry Smith MatNorm_SeqDense, 1987c0aa2d19SHong Zhang /* 20*/ MatAssemblyBegin_SeqDense, 1988c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1989905e6a2fSBarry Smith MatSetOption_SeqDense, 1990905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1991d519adbfSMatthew Knepley /* 24*/ MatZeroRows_SeqDense, 1992db4efbfdSBarry Smith 0, 1993db4efbfdSBarry Smith 0, 1994db4efbfdSBarry Smith 0, 1995db4efbfdSBarry Smith 0, 19964994cf47SJed Brown /* 29*/ MatSetUp_SeqDense, 1997273d9f13SBarry Smith 0, 1998905e6a2fSBarry Smith 0, 199973a71a0fSBarry Smith 0, 200073a71a0fSBarry Smith 0, 2001d519adbfSMatthew Knepley /* 34*/ MatDuplicate_SeqDense, 2002a5ae1ecdSBarry Smith 0, 2003a5ae1ecdSBarry Smith 0, 2004a5ae1ecdSBarry Smith 0, 2005a5ae1ecdSBarry Smith 0, 2006d519adbfSMatthew Knepley /* 39*/ MatAXPY_SeqDense, 2007a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 2008a5ae1ecdSBarry Smith 0, 20094b0e389bSBarry Smith MatGetValues_SeqDense, 2010a5ae1ecdSBarry Smith MatCopy_SeqDense, 2011d519adbfSMatthew Knepley /* 44*/ MatGetRowMax_SeqDense, 2012a5ae1ecdSBarry Smith MatScale_SeqDense, 2013a5ae1ecdSBarry Smith 0, 2014a5ae1ecdSBarry Smith 0, 2015a5ae1ecdSBarry Smith 0, 201673a71a0fSBarry Smith /* 49*/ MatSetRandom_SeqDense, 2017a5ae1ecdSBarry Smith 0, 2018a5ae1ecdSBarry Smith 0, 2019a5ae1ecdSBarry Smith 0, 2020a5ae1ecdSBarry Smith 0, 2021d519adbfSMatthew Knepley /* 54*/ 0, 2022a5ae1ecdSBarry Smith 0, 2023a5ae1ecdSBarry Smith 0, 2024a5ae1ecdSBarry Smith 0, 2025a5ae1ecdSBarry Smith 0, 2026d519adbfSMatthew Knepley /* 59*/ 0, 2027e03a110bSBarry Smith MatDestroy_SeqDense, 2028e03a110bSBarry Smith MatView_SeqDense, 2029357abbc8SBarry Smith 0, 203097304618SKris Buschelman 0, 2031d519adbfSMatthew Knepley /* 64*/ 0, 203297304618SKris Buschelman 0, 203397304618SKris Buschelman 0, 203497304618SKris Buschelman 0, 203597304618SKris Buschelman 0, 2036d519adbfSMatthew Knepley /* 69*/ MatGetRowMaxAbs_SeqDense, 203797304618SKris Buschelman 0, 203897304618SKris Buschelman 0, 203997304618SKris Buschelman 0, 204097304618SKris Buschelman 0, 2041d519adbfSMatthew Knepley /* 74*/ 0, 204297304618SKris Buschelman 0, 204397304618SKris Buschelman 0, 204497304618SKris Buschelman 0, 204597304618SKris Buschelman 0, 2046d519adbfSMatthew Knepley /* 79*/ 0, 204797304618SKris Buschelman 0, 204897304618SKris Buschelman 0, 204997304618SKris Buschelman 0, 20505bba2384SShri Abhyankar /* 83*/ MatLoad_SeqDense, 2051865e5f61SKris Buschelman 0, 20521cbb95d3SBarry Smith MatIsHermitian_SeqDense, 2053865e5f61SKris Buschelman 0, 2054865e5f61SKris Buschelman 0, 2055865e5f61SKris Buschelman 0, 2056d519adbfSMatthew Knepley /* 89*/ MatMatMult_SeqDense_SeqDense, 2057a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 2058a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 2059865e5f61SKris Buschelman 0, 2060865e5f61SKris Buschelman 0, 2061d519adbfSMatthew Knepley /* 94*/ 0, 20625df89d91SHong Zhang 0, 20635df89d91SHong Zhang 0, 20645df89d91SHong Zhang 0, 2065284134d9SBarry Smith 0, 2066d519adbfSMatthew Knepley /* 99*/ 0, 2067284134d9SBarry Smith 0, 2068284134d9SBarry Smith 0, 2069ba337c44SJed Brown MatConjugate_SeqDense, 2070f73d5cc4SBarry Smith 0, 2071ba337c44SJed Brown /*104*/ 0, 2072ba337c44SJed Brown MatRealPart_SeqDense, 2073ba337c44SJed Brown MatImaginaryPart_SeqDense, 2074985db425SBarry Smith 0, 2075985db425SBarry Smith 0, 207685e2c93fSHong Zhang /*109*/ MatMatSolve_SeqDense, 2077985db425SBarry Smith 0, 20788d0534beSBarry Smith MatGetRowMin_SeqDense, 2079aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 2080aabbc4fbSShri Abhyankar 0, 2081aabbc4fbSShri Abhyankar /*114*/ 0, 2082aabbc4fbSShri Abhyankar 0, 2083aabbc4fbSShri Abhyankar 0, 2084aabbc4fbSShri Abhyankar 0, 2085aabbc4fbSShri Abhyankar 0, 2086aabbc4fbSShri Abhyankar /*119*/ 0, 2087aabbc4fbSShri Abhyankar 0, 2088aabbc4fbSShri Abhyankar 0, 20890716a85fSBarry Smith 0, 20900716a85fSBarry Smith 0, 20910716a85fSBarry Smith /*124*/ 0, 20925df89d91SHong Zhang MatGetColumnNorms_SeqDense, 20935df89d91SHong Zhang 0, 20945df89d91SHong Zhang 0, 20955df89d91SHong Zhang 0, 20965df89d91SHong Zhang /*129*/ 0, 209775648e8dSHong Zhang MatTransposeMatMult_SeqDense_SeqDense, 209875648e8dSHong Zhang MatTransposeMatMultSymbolic_SeqDense_SeqDense, 209975648e8dSHong Zhang MatTransposeMatMultNumeric_SeqDense_SeqDense, 21003964eb88SJed Brown 0, 21013964eb88SJed Brown /*134*/ 0, 21023964eb88SJed Brown 0, 21033964eb88SJed Brown 0, 21043964eb88SJed Brown 0, 21053964eb88SJed Brown 0, 21063964eb88SJed Brown /*139*/ 0, 2107f9426fe0SMark Adams 0, 21083964eb88SJed Brown 0 2109985db425SBarry Smith }; 211090ace30eSBarry Smith 21114a2ae208SSatish Balay #undef __FUNCT__ 21124a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 21134b828684SBarry Smith /*@C 2114fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2115d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2116d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2117289bc588SBarry Smith 2118db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2119db81eaa0SLois Curfman McInnes 212020563c6bSBarry Smith Input Parameters: 2121db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 21220c775827SLois Curfman McInnes . m - number of rows 212318f449edSLois Curfman McInnes . n - number of columns 21240298fd71SBarry Smith - data - optional location of matrix data in column major order. Set data=NULL for PETSc 2125dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 212620563c6bSBarry Smith 212720563c6bSBarry Smith Output Parameter: 212844cd7ae7SLois Curfman McInnes . A - the matrix 212920563c6bSBarry Smith 2130b259b22eSLois Curfman McInnes Notes: 213118f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 213218f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 21330298fd71SBarry Smith set data=NULL. 213418f449edSLois Curfman McInnes 2135027ccd11SLois Curfman McInnes Level: intermediate 2136027ccd11SLois Curfman McInnes 2137dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2138d65003e9SLois Curfman McInnes 213969b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues() 214020563c6bSBarry Smith @*/ 21417087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2142289bc588SBarry Smith { 2143dfbe8321SBarry Smith PetscErrorCode ierr; 21443b2fbd54SBarry Smith 21453a40ed3dSBarry Smith PetscFunctionBegin; 2146f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2147f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2148273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2149273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2150273d9f13SBarry Smith PetscFunctionReturn(0); 2151273d9f13SBarry Smith } 2152273d9f13SBarry Smith 21534a2ae208SSatish Balay #undef __FUNCT__ 2154afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2155273d9f13SBarry Smith /*@C 2156273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2157273d9f13SBarry Smith 2158273d9f13SBarry Smith Collective on MPI_Comm 2159273d9f13SBarry Smith 2160273d9f13SBarry Smith Input Parameters: 2161273d9f13SBarry Smith + A - the matrix 21620298fd71SBarry Smith - data - the array (or NULL) 2163273d9f13SBarry Smith 2164273d9f13SBarry Smith Notes: 2165273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2166273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2167284134d9SBarry Smith need not call this routine. 2168273d9f13SBarry Smith 2169273d9f13SBarry Smith Level: intermediate 2170273d9f13SBarry Smith 2171273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2172273d9f13SBarry Smith 217369b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA() 2174867c911aSBarry Smith 2175273d9f13SBarry Smith @*/ 21767087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2177273d9f13SBarry Smith { 21784ac538c5SBarry Smith PetscErrorCode ierr; 2179a23d5eceSKris Buschelman 2180a23d5eceSKris Buschelman PetscFunctionBegin; 21814ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2182a23d5eceSKris Buschelman PetscFunctionReturn(0); 2183a23d5eceSKris Buschelman } 2184a23d5eceSKris Buschelman 2185a23d5eceSKris Buschelman #undef __FUNCT__ 2186afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 21877087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2188a23d5eceSKris Buschelman { 2189273d9f13SBarry Smith Mat_SeqDense *b; 2190dfbe8321SBarry Smith PetscErrorCode ierr; 2191273d9f13SBarry Smith 2192273d9f13SBarry Smith PetscFunctionBegin; 2193273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2194a868139aSShri Abhyankar 219534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 219634ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 219734ef9618SShri Abhyankar 2198273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 219986d161a7SShri Abhyankar b->Mmax = B->rmap->n; 220086d161a7SShri Abhyankar b->Nmax = B->cmap->n; 220186d161a7SShri Abhyankar if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 220286d161a7SShri Abhyankar 22039e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 22049e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2205*e92229d0SSatish Balay ierr = PetscCalloc1((size_t)b->lda*b->Nmax,&b->v);CHKERRQ(ierr); 22063bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 22072205254eSKarl Rupp 22089e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2209273d9f13SBarry Smith } else { /* user-allocated storage */ 22109e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2211273d9f13SBarry Smith b->v = data; 2212273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2213273d9f13SBarry Smith } 22140450473dSBarry Smith B->assembled = PETSC_TRUE; 2215273d9f13SBarry Smith PetscFunctionReturn(0); 2216273d9f13SBarry Smith } 2217273d9f13SBarry Smith 22181b807ce4Svictorle #undef __FUNCT__ 22191b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 22201b807ce4Svictorle /*@C 22211b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 22221b807ce4Svictorle 22231b807ce4Svictorle Input parameter: 22241b807ce4Svictorle + A - the matrix 22251b807ce4Svictorle - lda - the leading dimension 22261b807ce4Svictorle 22271b807ce4Svictorle Notes: 2228867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 22291b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 22301b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 22311b807ce4Svictorle 22321b807ce4Svictorle Level: intermediate 22331b807ce4Svictorle 22341b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 22351b807ce4Svictorle 2236284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2237867c911aSBarry Smith 22381b807ce4Svictorle @*/ 22397087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 22401b807ce4Svictorle { 22411b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 224221a2c019SBarry Smith 22431b807ce4Svictorle PetscFunctionBegin; 2244e32f2f54SBarry Smith if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n); 22451b807ce4Svictorle b->lda = lda; 224621a2c019SBarry Smith b->changelda = PETSC_FALSE; 224721a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 22481b807ce4Svictorle PetscFunctionReturn(0); 22491b807ce4Svictorle } 22501b807ce4Svictorle 22510bad9183SKris Buschelman /*MC 2252fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 22530bad9183SKris Buschelman 22540bad9183SKris Buschelman Options Database Keys: 22550bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 22560bad9183SKris Buschelman 22570bad9183SKris Buschelman Level: beginner 22580bad9183SKris Buschelman 225989665df3SBarry Smith .seealso: MatCreateSeqDense() 226089665df3SBarry Smith 22610bad9183SKris Buschelman M*/ 22620bad9183SKris Buschelman 22634a2ae208SSatish Balay #undef __FUNCT__ 22644a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 22658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat B) 2266273d9f13SBarry Smith { 2267273d9f13SBarry Smith Mat_SeqDense *b; 2268dfbe8321SBarry Smith PetscErrorCode ierr; 22697c334f02SBarry Smith PetscMPIInt size; 2270273d9f13SBarry Smith 2271273d9f13SBarry Smith PetscFunctionBegin; 2272ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr); 2273e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 227455659b69SBarry Smith 2275b00a9115SJed Brown ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 2276549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 227744cd7ae7SLois Curfman McInnes B->data = (void*)b; 227818f449edSLois Curfman McInnes 227944cd7ae7SLois Curfman McInnes b->pivots = 0; 2280273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2281273d9f13SBarry Smith b->v = 0; 228221a2c019SBarry Smith b->changelda = PETSC_FALSE; 22834e220ebcSLois Curfman McInnes 2284bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense);CHKERRQ(ierr); 2285bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr); 2286b24902e0SBarry Smith 2287bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr); 2288bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2289bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 2290bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 2291bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 2292bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 22933bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMult_seqaij_seqdense_C",MatTransposeMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 22943bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultSymbolic_seqaij_seqdense_C",MatTransposeMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 22953bf78175SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatTransposeMatMultNumeric_seqaij_seqdense_C",MatTransposeMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 229617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 22973a40ed3dSBarry Smith PetscFunctionReturn(0); 2298289bc588SBarry Smith } 2299