Actual source code: matregis.c
petsc-master 2020-08-25
2: #include <petsc/private/matimpl.h>
4: PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat);
5: PETSC_EXTERN PetscErrorCode MatCreate_MAIJ(Mat);
6: PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat);
7: PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat);
9: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
10: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
12: PETSC_EXTERN PetscErrorCode MatCreate_SeqBAIJ(Mat);
13: PETSC_EXTERN PetscErrorCode MatCreate_MPIBAIJ(Mat);
15: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat);
16: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat);
18: PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);
19: PETSC_EXTERN PetscErrorCode MatCreate_MPIDense(Mat);
20: #if defined(PETSC_HAVE_CUDA)
21: PETSC_EXTERN PetscErrorCode MatCreate_SeqDenseCUDA(Mat);
22: PETSC_EXTERN PetscErrorCode MatCreate_MPIDenseCUDA(Mat);
23: #endif
25: PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat);
26: PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat);
27: PETSC_EXTERN PetscErrorCode MatCreate_Composite(Mat);
29: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJPERM(Mat);
30: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJPERM(Mat);
32: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJSELL(Mat);
33: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJSELL(Mat);
35: #if defined(PETSC_HAVE_MKL_SPARSE)
36: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJMKL(Mat);
37: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJMKL(Mat);
38: #endif
40: #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
41: PETSC_EXTERN PetscErrorCode MatCreate_SeqBAIJMKL(Mat);
42: PETSC_EXTERN PetscErrorCode MatCreate_MPIBAIJMKL(Mat);
43: #endif
45: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCRL(Mat);
46: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJCRL(Mat);
48: PETSC_EXTERN PetscErrorCode MatCreate_Scatter(Mat);
49: PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat);
50: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat);
52: PETSC_EXTERN PetscErrorCode MatCreate_SeqSELL(Mat);
53: PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat);
55: #if defined(PETSC_HAVE_CUDA)
56: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJCUSPARSE(Mat);
57: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJCUSPARSE(Mat);
58: #endif
60: #if defined(PETSC_HAVE_VIENNACL)
61: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJViennaCL(Mat);
62: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJViennaCL(Mat);
63: #endif
65: #if defined(PETSC_HAVE_FFTW)
66: PETSC_EXTERN PetscErrorCode MatCreate_FFTW(Mat);
67: #endif
68: #if defined(PETSC_HAVE_ELEMENTAL)
69: PETSC_EXTERN PetscErrorCode MatCreate_Elemental(Mat);
70: #endif
71: #if defined(PETSC_HAVE_SCALAPACK)
72: PETSC_EXTERN PetscErrorCode MatCreate_ScaLAPACK(Mat);
73: #endif
75: PETSC_EXTERN PetscErrorCode MatCreate_Preallocator(Mat);
76: PETSC_EXTERN PetscErrorCode MatCreate_Dummy(Mat);
78: #if defined(PETSC_HAVE_HYPRE)
79: PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat);
80: #endif
82: PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat);
84: #if defined(PETSC_HAVE_HARA)
85: PETSC_EXTERN PetscErrorCode MatCreate_HARA(Mat);
86: #endif
88: /*@C
89: MatRegisterAll - Registers all of the matrix types in PETSc
91: Not Collective
93: Level: advanced
95: .seealso: MatRegister()
96: @*/
97: PetscErrorCode MatRegisterAll(void)
98: {
102: if (MatRegisterAllCalled) return(0);
103: MatRegisterAllCalled = PETSC_TRUE;
105: MatRegister(MATMFFD, MatCreate_MFFD);
107: MatRegister(MATMPIMAIJ, MatCreate_MAIJ);
108: MatRegister(MATSEQMAIJ, MatCreate_MAIJ);
109: MatRegister(MATMAIJ, MatCreate_MAIJ);
111: MatRegister(MATMPIKAIJ, MatCreate_KAIJ);
112: MatRegister(MATSEQKAIJ, MatCreate_KAIJ);
113: MatRegister(MATKAIJ, MatCreate_KAIJ);
115: MatRegister(MATIS, MatCreate_IS);
116: MatRegister(MATSHELL, MatCreate_Shell);
117: MatRegister(MATCOMPOSITE, MatCreate_Composite);
119: MatRegisterRootName(MATAIJ,MATSEQAIJ,MATMPIAIJ);
120: MatRegister(MATMPIAIJ, MatCreate_MPIAIJ);
121: MatRegister(MATSEQAIJ, MatCreate_SeqAIJ);
123: MatRegisterRootName(MATAIJPERM,MATSEQAIJPERM,MATMPIAIJPERM);
124: MatRegister(MATMPIAIJPERM, MatCreate_MPIAIJPERM);
125: MatRegister(MATSEQAIJPERM, MatCreate_SeqAIJPERM);
127: MatRegisterRootName(MATAIJSELL,MATSEQAIJSELL,MATMPIAIJSELL);
128: MatRegister(MATMPIAIJSELL, MatCreate_MPIAIJSELL);
129: MatRegister(MATSEQAIJSELL, MatCreate_SeqAIJSELL);
131: #if defined(PETSC_HAVE_MKL_SPARSE)
132: MatRegisterRootName(MATAIJMKL, MATSEQAIJMKL,MATMPIAIJMKL);
133: MatRegister(MATMPIAIJMKL, MatCreate_MPIAIJMKL);
134: MatRegister(MATSEQAIJMKL, MatCreate_SeqAIJMKL);
135: #endif
137: #if defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE)
138: MatRegisterRootName(MATBAIJMKL,MATSEQBAIJMKL,MATMPIBAIJMKL);
139: MatRegister(MATMPIBAIJMKL, MatCreate_MPIBAIJMKL);
140: MatRegister(MATSEQBAIJMKL, MatCreate_SeqBAIJMKL);
141: #endif
143: MatRegisterRootName(MATAIJCRL,MATSEQAIJCRL,MATMPIAIJCRL);
144: MatRegister(MATSEQAIJCRL, MatCreate_SeqAIJCRL);
145: MatRegister(MATMPIAIJCRL, MatCreate_MPIAIJCRL);
147: MatRegisterRootName(MATBAIJ,MATSEQBAIJ,MATMPIBAIJ);
148: MatRegister(MATMPIBAIJ, MatCreate_MPIBAIJ);
149: MatRegister(MATSEQBAIJ, MatCreate_SeqBAIJ);
151: MatRegisterRootName(MATSBAIJ,MATSEQSBAIJ,MATMPISBAIJ);
152: MatRegister(MATMPISBAIJ, MatCreate_MPISBAIJ);
153: MatRegister(MATSEQSBAIJ, MatCreate_SeqSBAIJ);
155: MatRegisterRootName(MATDENSE,MATSEQDENSE,MATMPIDENSE);
156: MatRegister(MATMPIDENSE, MatCreate_MPIDense);
157: MatRegister(MATSEQDENSE, MatCreate_SeqDense);
158: #if defined(PETSC_HAVE_CUDA)
159: MatRegisterRootName(MATDENSECUDA,MATSEQDENSECUDA,MATMPIDENSECUDA);
160: MatRegister(MATSEQDENSECUDA, MatCreate_SeqDenseCUDA);
161: MatRegister(MATMPIDENSECUDA, MatCreate_MPIDenseCUDA);
162: #endif
164: MatRegister(MATMPIADJ, MatCreate_MPIAdj);
165: MatRegister(MATSCATTER, MatCreate_Scatter);
166: MatRegister(MATBLOCKMAT, MatCreate_BlockMat);
167: MatRegister(MATNEST, MatCreate_Nest);
169: MatRegisterRootName(MATSELL,MATSEQSELL,MATMPISELL);
170: MatRegister(MATMPISELL, MatCreate_MPISELL);
171: MatRegister(MATSEQSELL, MatCreate_SeqSELL);
173: #if defined(PETSC_HAVE_CUDA)
174: MatRegisterRootName(MATAIJCUSPARSE,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE);
175: MatRegister(MATSEQAIJCUSPARSE, MatCreate_SeqAIJCUSPARSE);
176: MatRegister(MATMPIAIJCUSPARSE, MatCreate_MPIAIJCUSPARSE);
177: #endif
179: #if defined(PETSC_HAVE_VIENNACL)
180: MatRegisterRootName(MATAIJVIENNACL,MATSEQAIJVIENNACL,MATMPIAIJVIENNACL);
181: MatRegister(MATSEQAIJVIENNACL, MatCreate_SeqAIJViennaCL);
182: MatRegister(MATMPIAIJVIENNACL, MatCreate_MPIAIJViennaCL);
183: #endif
185: #if defined(PETSC_HAVE_FFTW)
186: MatRegister(MATFFTW, MatCreate_FFTW);
187: #endif
188: #if defined(PETSC_HAVE_ELEMENTAL)
189: MatRegister(MATELEMENTAL, MatCreate_Elemental);
190: #endif
191: #if defined(PETSC_HAVE_SCALAPACK)
192: MatRegister(MATSCALAPACK, MatCreate_ScaLAPACK);
193: #endif
195: MatRegister(MATPREALLOCATOR, MatCreate_Preallocator);
196: MatRegister(MATDUMMY, MatCreate_Dummy);
198: MatRegister(MATCONSTANTDIAGONAL,MatCreate_ConstantDiagonal);
200: #if defined(PETSC_HAVE_HYPRE)
201: MatRegister(MATHYPRE, MatCreate_HYPRE);
202: #endif
204: #if defined(PETSC_HAVE_HARA)
205: MatRegister(MATHARA, MatCreate_HARA);
206: #endif
207: return(0);
208: }