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-2011, Universitat 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> /*I "slepcst.h" I*/
26: typedef struct {
27: PetscBool setmat;
28: } ST_PRECOND;
32: PetscErrorCode STSetFromOptions_Precond(ST st)
33: {
35: PC pc;
36: const PCType pctype;
37: Mat P;
38: PetscBool t0,t1;
39: KSP ksp;
42: STGetKSP(st,&ksp);
43: KSPGetPC(ksp,&pc);
44: PetscObjectGetType((PetscObject)pc,&pctype);
45: STPrecondGetMatForPC(st,&P);
46: if (!pctype && st->A) {
47: if (P || st->shift_matrix == ST_MATMODE_SHELL) {
48: PCSetType(pc,PCJACOBI);
49: } else {
50: MatHasOperation(st->A,MATOP_DUPLICATE,&t0);
51: if (st->B) {
52: MatHasOperation(st->A,MATOP_AXPY,&t1);
53: } else {
54: t1 = PETSC_TRUE;
55: }
56: PCSetType(pc,(t0 && t1)?PCJACOBI:PCNONE);
57: }
58: }
59: return(0);
60: }
64: PetscErrorCode STSetUp_Precond(ST st)
65: {
66: Mat P;
67: PC pc;
68: PetscBool t0,setmat,destroyP=PETSC_FALSE,builtP;
72: /* if the user did not set the shift, use the target value */
73: if (!st->sigma_set) st->sigma = st->defsigma;
75: /* If pc is none and no matrix has to be set, exit */
76: STSetFromOptions_Precond(st);
77: if (!st->ksp) { STGetKSP(st,&st->ksp); }
78: KSPGetPC(st->ksp,&pc);
79: PetscTypeCompare((PetscObject)pc,PCNONE,&t0);
80: STPrecondGetKSPHasMat(st,&setmat);
81: if (t0 && !setmat) return(0);
83: /* Check if a user matrix is set */
84: STPrecondGetMatForPC(st,&P);
86: /* If not, create A - shift*B */
87: if (P) {
88: builtP = PETSC_FALSE;
89: destroyP = PETSC_TRUE;
90: PetscObjectReference((PetscObject)P);
91: } else {
92: builtP = PETSC_TRUE;
94: if (!(PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) && st->B) {
95: P = st->B;
96: destroyP = PETSC_FALSE;
97: } else if (st->sigma == 0.0) {
98: P = st->A;
99: destroyP = PETSC_FALSE;
100: } else if (PetscAbsScalar(st->sigma) < PETSC_MAX_REAL && st->shift_matrix != ST_MATMODE_SHELL) {
101: if (st->shift_matrix == ST_MATMODE_INPLACE) {
102: P = st->A;
103: destroyP = PETSC_FALSE;
104: } else {
105: MatDuplicate(st->A,MAT_COPY_VALUES,&P);
106: destroyP = PETSC_TRUE;
107: }
108: if (st->B) {
109: MatAXPY(P,-st->sigma,st->B,st->str);
110: } else {
111: MatShift(P,-st->sigma);
112: }
113: } else
114: builtP = PETSC_FALSE;
115: }
117: /* If P was not possible to obtain, set pc to PCNONE */
118: if (!P) {
119: PCSetType(pc,PCNONE);
121: /* If some matrix has to be set to ksp, a shell matrix is created */
122: if (setmat) {
123: STMatShellCreate(st,&P);
124: destroyP = PETSC_TRUE;
125: }
126: }
128: KSPSetOperators(st->ksp,setmat?P:PETSC_NULL,P,DIFFERENT_NONZERO_PATTERN);
130: if (destroyP) {
131: MatDestroy(&P);
132: } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
133: if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
134: if (st->B) {
135: MatAXPY(st->A,st->sigma,st->B,st->str);
136: } else {
137: MatShift(st->A,st->sigma);
138: }
139: }
140: }
141: return(0);
142: }
146: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
147: {
151: /* Nothing to be done if STSetUp has not been called yet */
152: if (!st->setupcalled) return(0);
153: st->sigma = newshift;
154: if (st->shift_matrix != ST_MATMODE_SHELL) {
155: STSetUp_Precond(st);
156: }
157: return(0);
158: }
163: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
164: {
166: PC pc;
167: PetscBool flag;
170: if (!st->ksp) { STGetKSP(st,&st->ksp); }
171: KSPGetPC(st->ksp,&pc);
172: PCGetOperatorsSet(pc,PETSC_NULL,&flag);
173: if (flag) {
174: PCGetOperators(pc,PETSC_NULL,mat,PETSC_NULL);
175: } else *mat = PETSC_NULL;
176: return(0);
177: }
182: /*@
183: STPrecondGetMatForPC - Returns the matrix previously set by STPrecondSetMatForPC().
185: Not Collective, but the Mat is shared by all processors that share the ST
187: Input Parameter:
188: . st - the spectral transformation context
190: Output Parameter:
191: . mat - the matrix that will be used in constructing the preconditioner or
192: PETSC_NULL if no matrix was set by STPrecondSetMatForPC().
194: Level: advanced
196: .seealso: STPrecondSetMatForPC()
197: @*/
198: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
199: {
205: PetscTryMethod(st,"STPrecondGetMatForPC_C",(ST,Mat*),(st,mat));
206: return(0);
207: }
212: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
213: {
214: PC pc;
215: Mat A;
216: PetscBool flag;
220: if (!st->ksp) { STGetKSP(st,&st->ksp); }
221: KSPGetPC(st->ksp,&pc);
222: /* Yes, all these lines are needed to safely set mat as the preconditioner
223: matrix in pc */
224: PCGetOperatorsSet(pc,&flag,PETSC_NULL);
225: if (flag) {
226: PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
227: PetscObjectReference((PetscObject)A);
228: } else
229: A = PETSC_NULL;
230: PetscObjectReference((PetscObject)mat);
231: PCSetOperators(pc,A,mat,DIFFERENT_NONZERO_PATTERN);
232: MatDestroy(&A);
233: MatDestroy(&mat);
234: STPrecondSetKSPHasMat(st,PETSC_TRUE);
235: return(0);
236: }
241: /*@
242: STPrecondSetMatForPC - Sets the matrix that must be used to build the preconditioner.
244: Logically Collective on ST and Mat
246: Input Parameter:
247: + st - the spectral transformation context
248: - mat - the matrix that will be used in constructing the preconditioner
250: Level: advanced
252: Notes:
253: This matrix will be passed to the KSP object (via KSPSetOperators) as
254: the matrix to be used when constructing the preconditioner.
255: If no matrix is set or mat is set to PETSC_NULL, A - sigma*B will
256: be used to build the preconditioner, being sigma the value set by STSetShift().
258: .seealso: STPrecondSetMatForPC(), STSetShift()
259: @*/
260: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
261: {
268: PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
269: return(0);
270: }
275: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscBool setmat)
276: {
277: ST_PRECOND *data = (ST_PRECOND*)st->data;
280: data->setmat = setmat;
281: return(0);
282: }
287: /*@
288: STPrecondSetKSPHasMat - Sets a flag indicating that during STSetUp the coefficient
289: matrix of the KSP linear system (A) must be set to be the same matrix as the
290: preconditioner (P).
292: Collective on ST
294: Input Parameter:
295: + st - the spectral transformation context
296: - setmat - the flag
298: Notes:
299: In most cases, the preconditioner matrix is used only in the PC object, but
300: in external solvers this matrix must be provided also as the A-matrix in
301: the KSP object.
303: Level: developer
305: .seealso: STPrecondGetKSPHasMat(), STSetShift()
306: @*/
307: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
308: {
314: PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
315: return(0);
316: }
321: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
322: {
323: ST_PRECOND *data = (ST_PRECOND*)st->data;
326: *setmat = data->setmat;
327: return(0);
328: }
333: /*@
334: STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
335: matrix of the KSP linear system (A) is set to be the same matrix as the
336: preconditioner (P).
338: Not Collective
340: Input Parameter:
341: . st - the spectral transformation context
343: Output Parameter:
344: . setmat - the flag
346: Level: developer
348: .seealso: STPrecondSetKSPHasMat(), STSetShift()
349: @*/
350: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
351: {
357: PetscTryMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
358: return(0);
359: }
363: PetscErrorCode STDestroy_Precond(ST st)
364: {
368: PetscFree(st->data);
369: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","",PETSC_NULL);
370: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","",PETSC_NULL);
371: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","",PETSC_NULL);
372: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","",PETSC_NULL);
373: return(0);
374: }
379: PetscErrorCode STCreate_Precond(ST st)
380: {
384: PetscNewLog(st,ST_PRECOND,&st->data);
385: st->ops->getbilinearform = STGetBilinearForm_Default;
386: st->ops->setup = STSetUp_Precond;
387: st->ops->setshift = STSetShift_Precond;
388: st->ops->destroy = STDestroy_Precond;
389: st->ops->setfromoptions = STSetFromOptions_Precond;
391: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","STPrecondGetMatForPC_Precond",STPrecondGetMatForPC_Precond);
392: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","STPrecondSetMatForPC_Precond",STPrecondSetMatForPC_Precond);
393: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","STPrecondGetKSPHasMat_Precond",STPrecondGetKSPHasMat_Precond);
394: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","STPrecondSetKSPHasMat_Precond",STPrecondSetKSPHasMat_Precond);
396: STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
397: return(0);
398: }