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-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepc-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: 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[0]) {
47: if (P || st->shift_matrix == ST_MATMODE_SHELL) {
48: PCSetType(pc,PCJACOBI);
49: } else {
50: MatHasOperation(st->A[0],MATOP_DUPLICATE,&t0);
51: if (st->nmat>1) {
52: MatHasOperation(st->A[0],MATOP_AXPY,&t1);
53: } else t1 = PETSC_TRUE;
54: PCSetType(pc,(t0 && t1)?PCJACOBI:PCNONE);
55: }
56: }
57: return(0);
58: }
62: PetscErrorCode STSetUp_Precond(ST st)
63: {
64: Mat P;
65: PC pc;
66: PetscBool t0,setmat,destroyP=PETSC_FALSE,builtP;
70: /* if the user did not set the shift, use the target value */
71: if (!st->sigma_set) st->sigma = st->defsigma;
73: /* If either pc is none and no matrix has to be set, or pc is shell , exit */
74: STSetFromOptions_Precond(st);
75: if (!st->ksp) { STGetKSP(st,&st->ksp); }
76: KSPGetPC(st->ksp,&pc);
77: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t0);
78: if (t0) return(0);
79: PetscObjectTypeCompare((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->nmat>1) {
95: P = st->A[1];
96: destroyP = PETSC_FALSE;
97: } else if (st->sigma == 0.0) {
98: P = st->A[0];
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[0];
103: destroyP = PETSC_FALSE;
104: } else {
105: MatDuplicate(st->A[0],MAT_COPY_VALUES,&P);
106: destroyP = PETSC_TRUE;
107: }
108: if (st->nmat>1) {
109: MatAXPY(P,-st->sigma,st->A[1],st->str);
110: } else {
111: MatShift(P,-st->sigma);
112: }
113: /* TODO: in case of ST_MATMODE_INPLACE should keep the Hermitian flag of st->A and restore at the end */
114: STMatSetHermitian(st,P);
115: } else builtP = PETSC_FALSE;
116: }
118: /* If P was not possible to obtain, set pc to PCNONE */
119: if (!P) {
120: PCSetType(pc,PCNONE);
122: /* If some matrix has to be set to ksp, a shell matrix is created */
123: if (setmat) {
124: STMatShellCreate(st,-st->sigma,0,NULL,&P);
125: STMatSetHermitian(st,P);
126: destroyP = PETSC_TRUE;
127: }
128: }
130: KSPSetOperators(st->ksp,setmat?P:NULL,P,DIFFERENT_NONZERO_PATTERN);
132: if (destroyP) {
133: MatDestroy(&P);
134: } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
135: if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
136: if (st->nmat>1) {
137: MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
138: } else {
139: MatShift(st->A[0],st->sigma);
140: }
141: }
142: }
143: return(0);
144: }
148: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
149: {
153: /* Nothing to be done if STSetUp has not been called yet */
154: if (!st->setupcalled) return(0);
155: st->sigma = newshift;
156: if (st->shift_matrix != ST_MATMODE_SHELL) {
157: STSetUp_Precond(st);
158: }
159: return(0);
160: }
164: static PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
165: {
167: PC pc;
168: PetscBool flag;
171: if (!st->ksp) { STGetKSP(st,&st->ksp); }
172: KSPGetPC(st->ksp,&pc);
173: PCGetOperatorsSet(pc,NULL,&flag);
174: if (flag) {
175: PCGetOperators(pc,NULL,mat,NULL);
176: } else *mat = NULL;
177: return(0);
178: }
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: 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: }
211: static PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
212: {
213: PC pc;
214: Mat A;
215: PetscBool flag;
219: if (!st->ksp) { STGetKSP(st,&st->ksp); }
220: KSPGetPC(st->ksp,&pc);
221: /* Yes, all these lines are needed to safely set mat as the preconditioner
222: matrix in pc */
223: PCGetOperatorsSet(pc,&flag,NULL);
224: if (flag) {
225: PCGetOperators(pc,&A,NULL,NULL);
226: PetscObjectReference((PetscObject)A);
227: } else A = NULL;
228: PetscObjectReference((PetscObject)mat);
229: PCSetOperators(pc,A,mat,DIFFERENT_NONZERO_PATTERN);
230: MatDestroy(&A);
231: MatDestroy(&mat);
232: STPrecondSetKSPHasMat(st,PETSC_TRUE);
233: return(0);
234: }
238: /*@
239: STPrecondSetMatForPC - Sets the matrix that must be used to build the preconditioner.
241: Logically Collective on ST and Mat
243: Input Parameter:
244: + st - the spectral transformation context
245: - mat - the matrix that will be used in constructing the preconditioner
247: Level: advanced
249: Notes:
250: This matrix will be passed to the KSP object (via KSPSetOperators) as
251: the matrix to be used when constructing the preconditioner.
252: If no matrix is set or mat is set to NULL, A - sigma*B will
253: be used to build the preconditioner, being sigma the value set by STSetShift().
255: .seealso: STPrecondSetMatForPC(), STSetShift()
256: @*/
257: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
258: {
265: PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
266: return(0);
267: }
271: static PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscBool setmat)
272: {
273: ST_PRECOND *data = (ST_PRECOND*)st->data;
276: data->setmat = setmat;
277: return(0);
278: }
282: /*@
283: STPrecondSetKSPHasMat - Sets a flag indicating that during STSetUp the coefficient
284: matrix of the KSP linear system (A) must be set to be the same matrix as the
285: preconditioner (P).
287: Collective on ST
289: Input Parameter:
290: + st - the spectral transformation context
291: - setmat - the flag
293: Notes:
294: In most cases, the preconditioner matrix is used only in the PC object, but
295: in external solvers this matrix must be provided also as the A-matrix in
296: the KSP object.
298: Level: developer
300: .seealso: STPrecondGetKSPHasMat(), STSetShift()
301: @*/
302: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
303: {
309: PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
310: return(0);
311: }
315: static PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
316: {
317: ST_PRECOND *data = (ST_PRECOND*)st->data;
320: *setmat = data->setmat;
321: return(0);
322: }
326: /*@
327: STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
328: matrix of the KSP linear system (A) is set to be the same matrix as the
329: preconditioner (P).
331: Not Collective
333: Input Parameter:
334: . st - the spectral transformation context
336: Output Parameter:
337: . setmat - the flag
339: Level: developer
341: .seealso: STPrecondSetKSPHasMat(), STSetShift()
342: @*/
343: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
344: {
350: PetscTryMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
351: return(0);
352: }
356: PetscErrorCode STDestroy_Precond(ST st)
357: {
361: PetscFree(st->data);
362: PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",NULL);
363: PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",NULL);
364: PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",NULL);
365: PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",NULL);
366: return(0);
367: }
371: PETSC_EXTERN PetscErrorCode STCreate_Precond(ST st)
372: {
376: PetscNewLog(st,ST_PRECOND,&st->data);
377: st->ops->getbilinearform = STGetBilinearForm_Default;
378: st->ops->setup = STSetUp_Precond;
379: st->ops->setshift = STSetShift_Precond;
380: st->ops->destroy = STDestroy_Precond;
381: st->ops->setfromoptions = STSetFromOptions_Precond;
383: PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",STPrecondGetMatForPC_Precond);
384: PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",STPrecondSetMatForPC_Precond);
385: PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",STPrecondGetKSPHasMat_Precond);
386: PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",STPrecondSetKSPHasMat_Precond);
388: STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
389: return(0);
390: }