Actual source code: dense.h

petsc-3.14.6 2021-03-30
Report Typos and Errors

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

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

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

 25:   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
 26:   Mat               cmat;      /* matrix representation of a given subset of columns */
 27:   Vec               cvec;      /* vector representation of a given column */
 28:   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
 29:   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
 30:   PetscInt          matinuse;  /* if cmat is in use (cbegin = matinuse-1) */
 31: } Mat_SeqDense;

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

 48: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat);
 49: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat);
 50: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat);
 51: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat);

 53: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);

 55: /* Used by SeqDenseCUDA */
 56: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
 57: PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
 58: PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
 59: PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
 60: PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
 61: PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
 62: PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
 63: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
 64: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
 65: PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
 66: PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
 67: PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
 68: PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
 69: PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
 70: PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
 71: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
 72: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
 73: PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);
 74: PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat,Vec,PetscInt);
 75: PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat,PetscScalar);
 76: PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat,PetscInt,Vec*);
 77: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat,PetscInt,Vec*);
 78: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
 79: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
 80: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
 81: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
 82: PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat,PetscInt,PetscInt,Mat*);
 83: PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat,Mat*);
 84: PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat,PetscInt);

 86: #if defined(PETSC_HAVE_CUDA)
 87: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
 88: PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
 89: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
 90: #endif

 92: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);

 94: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
 95: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

 97: PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer);
 98: PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer);
 99: PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
100: #endif