Actual source code: sbaijspooles.c
1: /*
2: Provides an interface to the Spooles serial sparse solver
3: */
5: #include src/mat/impls/aij/seq/spooles/spooles.h
9: PetscErrorCode MatDestroy_SeqSBAIJSpooles(Mat A)
10: {
12:
14: /* SeqSBAIJ_Spooles isn't really the matrix that USES spooles, */
15: /* rather it is a factory class for creating a symmetric matrix that can */
16: /* invoke Spooles' sequential cholesky solver. */
17: /* As a result, we don't have to clean up the stuff set by spooles */
18: /* as in MatDestroy_SeqAIJ_Spooles. */
19: MatConvert_Spooles_Base(A,MATSEQSBAIJ,&A);
20: (*A->ops->destroy)(A);
21: return(0);
22: }
26: PetscErrorCode MatAssemblyEnd_SeqSBAIJSpooles(Mat A,MatAssemblyType mode) \
27: {
29: PetscInt bs;
30: Mat_Spooles *lu=(Mat_Spooles *)(A->spptr);
33: (*lu->MatAssemblyEnd)(A,mode);
34: MatGetBlockSize(A,&bs);
35: if (bs > 1) SETERRQ1(PETSC_ERR_SUP,"Block size %D not supported by Spooles",bs);
36: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
37: A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJSpooles;
38: return(0);
39: }
41: /*
42: input:
43: F: numeric factor
44: output:
45: nneg, nzero, npos: matrix inertia
46: */
50: PetscErrorCode MatGetInertia_SeqSBAIJSpooles(Mat F,int *nneg,int *nzero,int *npos)
51: {
52: Mat_Spooles *lu = (Mat_Spooles*)F->spptr;
53: int neg,zero,pos;
56: FrontMtx_inertia(lu->frontmtx, &neg, &zero, &pos);
57: if(nneg) *nneg = neg;
58: if(nzero) *nzero = zero;
59: if(npos) *npos = pos;
60: return(0);
61: }
63: /* Note the Petsc r permutation is ignored */
66: PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJSpooles(Mat A,IS r,MatFactorInfo *info,Mat *F)
67: {
68: Mat B;
69: Mat_Spooles *lu;
71: int m=A->m,n=A->n;
74: /* Create the factorization matrix */
75: MatCreate(A->comm,m,n,m,n,&B);
76: MatSetType(B,A->type_name);
77: MatSeqSBAIJSetPreallocation(B,1,PETSC_NULL,PETSC_NULL);
79: B->ops->choleskyfactornumeric = MatFactorNumeric_SeqAIJSpooles;
80: B->ops->getinertia = MatGetInertia_SeqSBAIJSpooles;
81: B->factor = FACTOR_CHOLESKY;
83: lu = (Mat_Spooles *)(B->spptr);
84: lu->options.pivotingflag = SPOOLES_NO_PIVOTING;
85: lu->options.symflag = SPOOLES_SYMMETRIC; /* default */
86: lu->flg = DIFFERENT_NONZERO_PATTERN;
87: lu->options.useQR = PETSC_FALSE;
89: *F = B;
90: return(0);
91: }
96: PetscErrorCode MatConvert_SeqSBAIJ_SeqSBAIJSpooles(Mat A,const MatType type,Mat *newmat) {
97: /* This routine is only called to convert a MATSEQSBAIJ matrix */
98: /* to a MATSEQSBAIJSPOOLES matrix, so we will ignore 'MatType type'. */
100: Mat B=*newmat;
101: Mat_Spooles *lu;
104: if (B != A) {
105: /* This routine is inherited, so we know the type is correct. */
106: MatDuplicate(A,MAT_COPY_VALUES,&B);
107: }
109: PetscNew(Mat_Spooles,&lu);
110: B->spptr = (void*)lu;
112: lu->basetype = MATSEQSBAIJ;
113: lu->CleanUpSpooles = PETSC_FALSE;
114: lu->MatDuplicate = A->ops->duplicate;
115: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
116: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
117: lu->MatView = A->ops->view;
118: lu->MatAssemblyEnd = A->ops->assemblyend;
119: lu->MatDestroy = A->ops->destroy;
121: B->ops->duplicate = MatDuplicate_Spooles;
122: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJSpooles;
123: B->ops->assemblyend = MatAssemblyEnd_SeqSBAIJSpooles;
124: B->ops->destroy = MatDestroy_SeqSBAIJSpooles;
125: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaijspooles_seqsbaij_C",
126: "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
127: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqsbaijspooles_C",
128: "MatConvert_SeqSBAIJ_SeqSBAIJSpooles",MatConvert_SeqSBAIJ_SeqSBAIJSpooles);
129: PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJSPOOLES);
130: *newmat = B;
131: return(0);
132: }
135: /*MC
136: MATSEQSBAIJSPOOLES - MATSEQSBAIJSPOOLES = "seqsbaijspooles" - A matrix type providing direct solvers (Cholesky) for sequential symmetric
137: matrices via the external package Spooles.
139: If Spooles is installed (see the manual for
140: instructions on how to declare the existence of external packages),
141: a matrix type can be constructed which invokes Spooles solvers.
142: After calling MatCreate(...,A), simply call MatSetType(A,MATSEQSBAIJSPOOLES).
143: This matrix type is only supported for double precision real.
145: This matrix inherits from MATSEQSBAIJ. As a result, MatSeqSBAIJSetPreallocation is
146: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
147: the MATSEQSBAIJ type without data copy.
149: Options Database Keys:
150: + -mat_type seqsbaijspooles - sets the matrix type to seqsbaijspooles during calls to MatSetFromOptions()
151: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
152: . -mat_spooles_seed <seed> - random number seed used for ordering
153: . -mat_spooles_msglvl <msglvl> - message output level
154: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
155: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
156: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
157: . -mat_spooles_maxsize <n> - maximum size of a supernode
158: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
159: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
160: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
161: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
162: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
163: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
164: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
166: Level: beginner
168: .seealso: MATMPISBAIJSPOOLES, MATSEQAIJSPOOLES, MATMPIAIJSPOOLES, PCCHOLESKY
169: M*/
174: PetscErrorCode MatCreate_SeqSBAIJSpooles(Mat A)
175: {
179: /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ */
180: /* and SeqSBAIJSpooles types */
181: PetscObjectChangeTypeName((PetscObject)A,MATSEQSBAIJSPOOLES);
182: MatSetType(A,MATSEQSBAIJ);
183: MatConvert_SeqSBAIJ_SeqSBAIJSpooles(A,MATSEQSBAIJSPOOLES,&A);
184: return(0);
185: }
188: /*MC
189: MATSBAIJSPOOLES - MATSBAIJSPOOLES = "sbaijspooles" - A matrix type providing direct solvers (Cholesky) for sequential and parallel symmetric matrices via the external package Spooles.
191: If Spooles is installed (see the manual for
192: instructions on how to declare the existence of external packages),
193: a matrix type can be constructed which invokes Spooles solvers.
194: After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJSPOOLES).
195: This matrix type is only supported for double precision real.
197: This matrix inherits from MATSBAIJ. As a result, MatSeqSBAIJSetPreallocation and MatMPISBAIJSetPreallocation are
198: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
199: the MATSBAIJ type without data copy.
201: Options Database Keys:
202: + -mat_type sbaijspooles - sets the matrix type to sbaijspooles during calls to MatSetFromOptions()
203: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
204: . -mat_spooles_seed <seed> - random number seed used for ordering
205: . -mat_spooles_msglvl <msglvl> - message output level
206: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
207: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
208: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
209: . -mat_spooles_maxsize <n> - maximum size of a supernode
210: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
211: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
212: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
213: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
214: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
215: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
216: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
218: Level: beginner
220: .seealso: MATMPISBAIJSPOOLES, MATSEQAIJSPOOLES, MATMPIAIJSPOOLES, PCCHOLESKY
221: M*/
226: PetscErrorCode MatCreate_SBAIJSpooles(Mat A)
227: {
229: int size;
232: /* Change type name before calling MatSetType to force proper construction of SeqSBAIJSpooles or MPISBAIJSpooles */
233: PetscObjectChangeTypeName((PetscObject)A,MATSBAIJSPOOLES);
234: MPI_Comm_size(A->comm,&size);
235: if (size == 1) {
236: MatSetType(A,MATSEQSBAIJ);
237: MatConvert_SeqSBAIJ_SeqSBAIJSpooles(A,MATSEQSBAIJSPOOLES,&A);
238: } else {
239: MatSetType(A,MATMPISBAIJ);
240: MatConvert_MPISBAIJ_MPISBAIJSpooles(A,MATMPISBAIJSPOOLES,&A);
241: }
242:
243: return(0);
244: }