Actual source code: sbaij.h

petsc-3.10.5 2019-03-28
Report Typos and Errors

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

  7: /*
  8:   MATSEQSBAIJ format - Block compressed row storage. The i[] and j[]
  9:   arrays start at 0.
 10: */

 12: typedef struct {
 13:   SEQAIJHEADER(MatScalar);
 14:   SEQBAIJHEADER;
 15:   PetscInt         *inew;        /* pointer to beginning of each row of reordered matrix */
 16:   PetscInt         *jnew;        /* column values: jnew + i[k] is start of row k */
 17:   MatScalar        *anew;        /* nonzero diagonal and superdiagonal elements of reordered matrix */
 18:   PetscScalar      *solves_work; /* work space used in MatSolves */
 19:   PetscInt         solves_work_n; /* size of solves_work */
 20:   PetscInt         *a2anew;        /* map used for symm permutation */
 21:   PetscBool        permute;        /* if true, a non-trivial permutation is used for factorization */
 22:   PetscBool        ignore_ltriangular; /* if true, ignore the lower triangular values inserted by users */
 23:   PetscBool        getrow_utriangular; /* if true, MatGetRow_SeqSBAIJ() is enabled to get the upper part of the row */
 24:   Mat_SeqAIJ_Inode inode;
 25:   unsigned short   *jshort;
 26:   PetscBool        free_jshort;
 27: } Mat_SeqSBAIJ;

 29: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
 30: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
 31: PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,const MatFactorInfo*);
 32: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
 33: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
 34: PETSC_INTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
 35: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
 36: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
 37: PETSC_INTERN PetscErrorCode MatSeqSBAIJZeroOps_Private(Mat);
 38: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqSBAIJ(Mat,IS,IS,MatReuse,Mat*);
 39: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
 40: PETSC_INTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar);
 41: PETSC_INTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal*);
 42: PETSC_INTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscBool*);
 43: PETSC_INTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
 44: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
 45: PETSC_INTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo*);
 46: PETSC_INTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
 47: PETSC_INTERN PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat,Vec,PetscInt[]);
 48: PETSC_INTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
 49: PETSC_INTERN PetscErrorCode MatDestroy_SeqSBAIJ(Mat);
 50: PETSC_INTERN PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);

 52: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 53: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Mat,const MatFactorInfo*);
 54: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
 55: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);

 57: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
 58: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
 59: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
 60: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);

 62: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
 63: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
 64: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);
 65: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);

 67: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 68: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
 69: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
 70: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);

 72: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 73: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
 74: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
 75: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);

 77: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 78: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
 79: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
 80: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);

 82: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 83: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
 84: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
 85: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);

 87: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 88: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
 89: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
 90: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);

 92: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 93: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
 94: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
 95: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);

 97: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
 98: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
 99: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
100: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
101: PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
102: PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);

104: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat,const MatFactorInfo*);
105: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat,Mat,const MatFactorInfo*);
106: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat,const MatFactorInfo*);
107: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat,const MatFactorInfo*);
108: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat,const MatFactorInfo*);
109: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat,const MatFactorInfo*);
110: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat,const MatFactorInfo*);
111: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat,const MatFactorInfo*);

113: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
114: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
115: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
116: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_inplace(Mat,Vec,Vec);
117: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_inplace(Mat,Vec,Vec);
118: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_inplace(Mat,Vec,Vec);
119: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_inplace(Mat,Vec,Vec);
120: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_inplace(Mat,Vec,Vec);
121: PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_inplace(Mat,Vec,Vec);

123: PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1_inplace(Mat,Vecs,Vecs);
124: PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);

126: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
127: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
128: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
129: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
130: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
131: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
132: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
133: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);

135: PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat,Vec,Vec);

137: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
138: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
139: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
140: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
141: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
142: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
143: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
144: PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);

146: PETSC_INTERN PetscErrorCode MatSOR_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
147: PETSC_INTERN PetscErrorCode MatLoad_SeqSBAIJ(Mat,PetscViewer);
148: PETSC_INTERN PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool);

150: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat,Mat,PetscInt*);

152: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
153: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
154: /* required by mpisbaij.c */
155: PETSC_INTERN PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
156: PETSC_INTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
157: PETSC_INTERN PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
158: PETSC_INTERN PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
159: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
160: PETSC_INTERN PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*,Vec,Vec);

162: #endif