Actual source code: dense.h
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: PetscScalar *tau; /* scalar factors of QR factorization */
21: Vec qrrhs; /* RHS for solving with QR (solution vector can't hold copy of RHS) */
22: PetscBLASInt lda; /* Lapack leading dimension of data */
23: PetscBLASInt rank; /* numerical rank (of a QR factorized matrix) */
24: PetscBool user_alloc; /* true if the user provided the dense data */
25: PetscBool unplaced_user_alloc;
26: Mat ptapwork; /* workspace (SeqDense matrix) for PtAP */
28: /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
29: Mat cmat; /* matrix representation of a given subset of columns */
30: Vec cvec; /* vector representation of a given column */
31: const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
32: PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */
33: PetscInt matinuse; /* if cmat is in use (cbegin = matinuse-1) */
34: } Mat_SeqDense;
36: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
37: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
38: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
39: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
40: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
41: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
42: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
43: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
44: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
45: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat);
46: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
47: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat);
48: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*);
49: PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat);
51: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat);
52: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat);
53: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat);
54: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat);
56: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);
58: /* Used by SeqDenseCUDA */
59: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
60: PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
61: PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
62: PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
63: PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
64: PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
65: PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
66: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
67: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
68: PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
69: PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
70: PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
71: PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
72: PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
73: PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
74: PETSC_INTERN PetscErrorCode MatQRFactor_SeqDense(Mat,IS,const MatFactorInfo*);
75: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
76: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
77: PETSC_INTERN PetscErrorCode MatQRFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
78: PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);
79: PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat,Vec,PetscInt);
80: PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat,PetscScalar);
81: PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat,PetscInt,Vec*);
82: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat,PetscInt,Vec*);
83: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
84: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
85: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
86: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
87: PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat,PetscInt,PetscInt,Mat*);
88: PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat,Mat*);
89: PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat,PetscInt);
90: PETSC_INTERN PetscErrorCode MatCopy_SeqDense(Mat,Mat,MatStructure);
91: PETSC_INTERN PetscErrorCode MatZeroEntries_SeqDense(Mat);
93: #if defined(PETSC_HAVE_CUDA)
94: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
95: PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
96: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
97: #endif
99: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
101: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
102: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
104: PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer);
105: PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer);
106: PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
107: #endif