Actual source code: dense.h

petsc-3.13.6 2020-09-29
Report Typos and Errors

  4:  #include <petsc/private/matimpl.h>
  5:  #include <../src/mat/impls/aij/seq/aij.h>

  7: /*
  8:   MATSEQDENSE format - conventional dense Fortran storage (by columns)
  9: */

 11: typedef struct {
 12:   PetscScalar  *v;                /* matrix elements */
 13:   PetscScalar  *unplacedarray;    /* if one called MatDensePlaceArray(), this is where it stashed the original */
 14:   PetscBool    roworiented;       /* if true, row oriented input (default) */
 15:   PetscInt     pad;               /* padding */
 16:   PetscBLASInt *pivots;           /* pivots in LU factorization */
 17:   PetscBLASInt lfwork;            /* length of work array in factorization */
 18:   PetscScalar  *fwork;            /* work array in factorization */
 19:   PetscBLASInt lda;               /* Lapack leading dimension of data */
 20:   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
 21:   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
 22:   PetscBool    user_alloc;        /* true if the user provided the dense data */
 23:   PetscBool    unplaced_user_alloc;
 24:   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */

 26:   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMultxxx_SeqAIJ_SeqDense */
 27: } Mat_SeqDense;

 29: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
 30: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
 31: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
 32: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
 33: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
 34: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
 35: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
 36: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
 37: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
 38: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat);
 39: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
 40: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat);
 41: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*);
 42: PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat);

 44: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat);
 45: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat);
 46: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat);
 47: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat);

 49: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);

 51: /* Used by SeqDenseCUDA */
 52: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
 53: PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
 54: PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
 55: PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
 56: PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
 57: PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
 58: PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
 59: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
 60: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
 61: PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
 62: PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
 63: PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
 64: PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
 65: PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
 66: PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
 67: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
 68: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
 69: PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);

 71: #if defined(PETSC_HAVE_CUDA)
 72: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
 73: PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
 74: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
 75: #endif

 77: PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
 78: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);

 80: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
 81: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

 83: PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer);
 84: PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer);
 85: PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
 86: #endif