Actual source code: precond.c
1: /*
2: Implements the ST class for preconditioned eigenvalue methods.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include private/stimpl.h
26: PetscErrorCode STDestroy_Precond(ST st);
27: PetscErrorCode STSetFromOptions_Precond(ST st);
29: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat);
30: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat);
31: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscTruth setmat);
32: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscTruth *setmat);
35: typedef struct {
36: PetscTruth setmat;
37: } ST_PRECOND;
42: PetscErrorCode SLEPcNotImplemented_Precond(ST st, Vec x, Vec y) {
43: SETERRQ(1, "STPrecond does not support some operation. Please, refer to the SLEPc Manual for more information.");
44: }
48: PetscErrorCode STSetFromOptions_Precond(ST st)
49: {
51: PC pc;
52: const PCType pctype;
53: Mat P;
54: PetscTruth t0, t1;
58: KSPGetPC(st->ksp, &pc);
59: PetscObjectGetType((PetscObject)pc, &pctype);
60: STPrecondGetMatForPC(st, &P);
61: if (!pctype && st->A) {
62: if (P || st->shift_matrix == ST_MATMODE_SHELL) {
63: PCSetType(pc, PCJACOBI);
64: } else {
65: MatHasOperation(st->A, MATOP_DUPLICATE, &t0);
66: if (st->B) {
67: MatHasOperation(st->A, MATOP_AXPY, &t1);
68: } else {
69: t1 = PETSC_TRUE;
70: }
71: PCSetType(pc, (t0 && t1)? PCJACOBI:PCNONE);
72: }
73: }
75: return(0);
76: }
81: PetscErrorCode STSetUp_Precond(ST st)
82: {
83: Mat P;
84: PC pc;
85: PetscTruth t0, setmat, destroyP=PETSC_FALSE, builtP;
90: /* if the user did not set the shift, use the target value */
91: if (!st->sigma_set) st->sigma = st->defsigma;
93: /* If pc is none and any matrix has to be set, exit */
94: STSetFromOptions_Precond(st);
95: KSPGetPC(st->ksp, &pc);
96: PetscTypeCompare((PetscObject)pc, PCNONE, &t0);
97: STPrecondGetKSPHasMat(st, &setmat);
98: if (t0 && !setmat) return(0);
100: /* Check if a user matrix is set */
101: STPrecondGetMatForPC(st, &P);
103: /* If not, create A - shift*B */
104: if (P) {
105: builtP = PETSC_FALSE;
106: destroyP = PETSC_TRUE;
107: PetscObjectReference((PetscObject)P);
108: } else {
109: builtP = PETSC_TRUE;
111: if (st->shift_matrix == ST_MATMODE_SHELL) {
112: STMatShellCreate(st,&P);
113: //TODO: set the apply and apply transpose to st->mat
114: destroyP = PETSC_TRUE;
115: } else if (!(PetscAbsScalar(st->sigma) < PETSC_MAX) && st->B) {
116: P = st->B;
117: destroyP = PETSC_FALSE;
118: } else if (st->sigma == 0.0) {
119: P = st->A;
120: destroyP = PETSC_FALSE;
121: } else if (PetscAbsScalar(st->sigma) < PETSC_MAX) {
122: if (st->shift_matrix == ST_MATMODE_INPLACE) {
123: P = st->A;
124: destroyP = PETSC_FALSE;
125: } else {
126: MatDuplicate(st->A, MAT_COPY_VALUES, &P);
127: destroyP = PETSC_TRUE;
128: }
129: if (st->B) {
130: MatAXPY(P, -st->sigma, st->B, st->str);
131: } else {
132: MatShift(P, -st->sigma);
133: }
134: } else
135: builtP = PETSC_FALSE;
136: }
138: /* If P was not possible to obtain, set pc to PCNONE */
139: if (!P) {
140: PCSetType(pc, PCNONE);
142: /* If some matrix has to be set to ksp, set ksp to KSPPREONLY */
143: if (setmat) {
144: STMatShellCreate(st, &P);
145: destroyP = PETSC_TRUE;
146: KSPSetType(st->ksp, KSPPREONLY);
147: }
148: }
150: KSPSetOperators(st->ksp,setmat?P:PETSC_NULL,P,DIFFERENT_NONZERO_PATTERN);
152: if (destroyP) {
153: MatDestroy(P);
154: } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
155: if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX) {
156: if (st->B) {
157: MatAXPY(st->A,st->sigma,st->B,st->str);
158: } else {
159: MatShift(st->A,st->sigma);
160: }
161: }
162: }
164: return(0);
165: }
169: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
170: {
175: /* Nothing to be done if STSetUp has not been called yet */
176: if (!st->setupcalled) return(0);
177:
178: st->sigma = newshift;
179: if (st->shift_matrix != ST_MATMODE_SHELL) {
180: STSetUp_Precond(st);
181: }
183: return(0);
184: }
189: PetscErrorCode STCreate_Precond(ST st)
190: {
192: ST_PRECOND *data;
196: PetscNew(ST_PRECOND, &data);
197: st->data = data;
199: st->ops->apply = SLEPcNotImplemented_Precond;
200: st->ops->getbilinearform = STGetBilinearForm_Default;
201: st->ops->applytrans = SLEPcNotImplemented_Precond;
202: st->ops->postsolve = PETSC_NULL;
203: st->ops->backtr = PETSC_NULL;
204: st->ops->setup = STSetUp_Precond;
205: st->ops->setshift = STSetShift_Precond;
206: st->ops->view = STView_Default;
207: st->ops->destroy = STDestroy_Precond;
208: st->ops->setfromoptions = STSetFromOptions_Precond;
209:
210: st->checknullspace = PETSC_NULL;
212: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","STPrecondGetMatForPC_Precond",STPrecondGetMatForPC_Precond);
213: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","STPrecondSetMatForPC_Precond",STPrecondSetMatForPC_Precond);
214: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","STPrecondGetKSPHasMat_Precond",STPrecondGetKSPHasMat_Precond);
215: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","STPrecondSetKSPHasMat_Precond",STPrecondSetKSPHasMat_Precond);
217: STPrecondSetKSPHasMat_Precond(st, PETSC_TRUE);
218: KSPSetType(st->ksp, KSPPREONLY);
220: return(0);
221: }
226: PetscErrorCode STDestroy_Precond(ST st)
227: {
232: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","",PETSC_NULL);
233: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","",PETSC_NULL);
234: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","",PETSC_NULL);
235: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","",PETSC_NULL);
237: PetscFree(st->data);
239: return(0);
240: }
244: /*@
245: STPrecondGetMatForPC - Gets the matrix previously set by STPrecondSetMatForPC.
246: This matrix will be passed as parameter in the KSPSetOperator function as
247: the matrix to be used in constructing the preconditioner.
249: Collective on ST
251: Input Parameter:
252: . st - the spectral transformation context
254: Output Parameter:
255: . mat - the matrix that will be used in constructing the preconditioner or
256: PETSC_NULL if any previous matrix was set by STPrecondSetMatForPC.
258: Level: advanced
260: .seealso: STPrecondSetMatForPC(), KSPSetOperator()
261: @*/
262: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
263: {
264: PetscErrorCode ierr, (*f)(ST,Mat*);
268: PetscObjectQueryFunction((PetscObject)st,"STPrecondGetMatForPC_C",(void (**)())&f);
269: if (f) {
270: (*f)(st,mat);
271: }
272: return(0);
273: }
276: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
277: {
279: PC pc;
280: PetscTruth flag;
285: KSPGetPC(st->ksp, &pc);
286: PCGetOperatorsSet(pc, PETSC_NULL, &flag);
287: if (flag) {
288: PCGetOperators(pc, PETSC_NULL, mat, PETSC_NULL);
289: } else
290: *mat = PETSC_NULL;
292: return(0);
293: }
298: /*@
299: STPrecondSetMatForPC - Sets the matrix that will be passed as parameter in
300: the KSPSetOperators function as the matrix to be used in constructing the
301: preconditioner. If any matrix is set or mat is PETSC_NULL, A - sigma*B will
302: be used, being sigma the value set by STSetShift
304: Collective on ST
306: Input Parameter:
307: + st - the spectral transformation context
308: - mat - the matrix that will be used in constructing the preconditioner
310: Level: advanced
312: .seealso: STPrecondSetMatForPC(), KSPSetOperators()
313: @*/
314: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
315: {
316: PetscErrorCode ierr, (*f)(ST,Mat);
321: PetscObjectQueryFunction((PetscObject)st,"STPrecondSetMatForPC_C",(void (**)())&f);
322: if (f) {
323: (*f)(st,mat);
324: }
325: return(0);
326: }
329: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
330: {
331: PC pc;
332: Mat A;
333: PetscTruth flag;
340: KSPGetPC(st->ksp, &pc);
341: /* Yes, all these lines are needed to safely set mat as the preconditioner
342: matrix in pc */
343: PCGetOperatorsSet(pc, &flag, PETSC_NULL);
344: if (flag) {
345: PCGetOperators(pc, &A, PETSC_NULL, PETSC_NULL);
346: PetscObjectReference((PetscObject)A);
347: } else
348: A = PETSC_NULL;
349: PetscObjectReference((PetscObject)mat);
350: PCSetOperators(pc, A, mat, DIFFERENT_NONZERO_PATTERN);
351: if (A) { MatDestroy(A); }
352: MatDestroy(mat);
354: return(0);
355: }
361: /*@
362: STPrecondSetKSPHasMat - Sets if during the STSetUp the KSP matrix associated
363: to the linear system is set with the matrix for building the preconditioner.
365: Collective on ST
367: Input Parameter:
368: + st - the spectral transformation context
369: - setmat - if true, the KSP matrix associated to linear system is set with
370: the matrix for building the preconditioner
372: Level: developer
374: .seealso: STPrecondGetKSPHasMat(), TSetShift(), KSPSetOperators()
375: @*/
376: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscTruth setmat)
377: {
378: PetscErrorCode ierr, (*f)(ST,PetscTruth);
382: PetscObjectQueryFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",(void (**)())&f);
383: if (f) {
384: (*f)(st,setmat);
385: }
386: return(0);
387: }
391: /*@
392: STPrecondGetKSPHasMat - Gets if during the STSetUp the KSP matrix associated
393: to the linear system is set with the matrix for building the preconditioner.
395: Collective on ST
397: Input Parameter:
398: . st - the spectral transformation context
400: Output Parameter:
401: . setmat - if true, the KSP matrix associated to linear system is set with
402: the matrix for building the preconditioner
404: Level: developer
406: .seealso: STPrecondSetKSPHasMat(), STSetShift(), KSPSetOperators()
407: @*/
408: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscTruth *setmat)
409: {
410: PetscErrorCode ierr, (*f)(ST,PetscTruth*);
414: PetscObjectQueryFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",(void (**)())&f);
415: if (f) {
416: (*f)(st,setmat);
417: }
418: return(0);
419: }
422: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscTruth setmat)
423: {
424: ST_PRECOND *data = (ST_PRECOND*)st->data;
430: data->setmat = setmat;
432: return(0);
433: }
437: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscTruth *setmat)
438: {
439: ST_PRECOND *data = (ST_PRECOND*)st->data;
445: *setmat = data->setmat;
447: return(0);
448: }