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-2012, 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 <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: 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: 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->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: /* 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
116: builtP = PETSC_FALSE;
117: }
119: /* If P was not possible to obtain, set pc to PCNONE */
120: if (!P) {
121: PCSetType(pc,PCNONE);
123: /* If some matrix has to be set to ksp, a shell matrix is created */
124: if (setmat) {
125: STMatShellCreate(st,&P);
126: STMatSetHermitian(st,P);
127: destroyP = PETSC_TRUE;
128: }
129: }
131: KSPSetOperators(st->ksp,setmat?P:PETSC_NULL,P,DIFFERENT_NONZERO_PATTERN);
133: if (destroyP) {
134: MatDestroy(&P);
135: } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
136: if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
137: if (st->B) {
138: MatAXPY(st->A,st->sigma,st->B,st->str);
139: } else {
140: MatShift(st->A,st->sigma);
141: }
142: }
143: }
144: return(0);
145: }
149: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
150: {
154: /* Nothing to be done if STSetUp has not been called yet */
155: if (!st->setupcalled) return(0);
156: st->sigma = newshift;
157: if (st->shift_matrix != ST_MATMODE_SHELL) {
158: STSetUp_Precond(st);
159: }
160: return(0);
161: }
163: EXTERN_C_BEGIN
166: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
167: {
169: PC pc;
170: PetscBool flag;
173: if (!st->ksp) { STGetKSP(st,&st->ksp); }
174: KSPGetPC(st->ksp,&pc);
175: PCGetOperatorsSet(pc,PETSC_NULL,&flag);
176: if (flag) {
177: PCGetOperators(pc,PETSC_NULL,mat,PETSC_NULL);
178: } else *mat = PETSC_NULL;
179: return(0);
180: }
181: EXTERN_C_END
185: /*@
186: STPrecondGetMatForPC - Returns the matrix previously set by STPrecondSetMatForPC().
188: Not Collective, but the Mat is shared by all processors that share the ST
190: Input Parameter:
191: . st - the spectral transformation context
193: Output Parameter:
194: . mat - the matrix that will be used in constructing the preconditioner or
195: PETSC_NULL if no matrix was set by STPrecondSetMatForPC().
197: Level: advanced
199: .seealso: STPrecondSetMatForPC()
200: @*/
201: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
202: {
208: PetscTryMethod(st,"STPrecondGetMatForPC_C",(ST,Mat*),(st,mat));
209: return(0);
210: }
212: EXTERN_C_BEGIN
215: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
216: {
217: PC pc;
218: Mat A;
219: PetscBool flag;
223: if (!st->ksp) { STGetKSP(st,&st->ksp); }
224: KSPGetPC(st->ksp,&pc);
225: /* Yes, all these lines are needed to safely set mat as the preconditioner
226: matrix in pc */
227: PCGetOperatorsSet(pc,&flag,PETSC_NULL);
228: if (flag) {
229: PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
230: PetscObjectReference((PetscObject)A);
231: } else
232: A = PETSC_NULL;
233: PetscObjectReference((PetscObject)mat);
234: PCSetOperators(pc,A,mat,DIFFERENT_NONZERO_PATTERN);
235: MatDestroy(&A);
236: MatDestroy(&mat);
237: STPrecondSetKSPHasMat(st,PETSC_TRUE);
238: return(0);
239: }
240: EXTERN_C_END
244: /*@
245: STPrecondSetMatForPC - Sets the matrix that must be used to build the preconditioner.
247: Logically Collective on ST and Mat
249: Input Parameter:
250: + st - the spectral transformation context
251: - mat - the matrix that will be used in constructing the preconditioner
253: Level: advanced
255: Notes:
256: This matrix will be passed to the KSP object (via KSPSetOperators) as
257: the matrix to be used when constructing the preconditioner.
258: If no matrix is set or mat is set to PETSC_NULL, A - sigma*B will
259: be used to build the preconditioner, being sigma the value set by STSetShift().
261: .seealso: STPrecondSetMatForPC(), STSetShift()
262: @*/
263: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
264: {
271: PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
272: return(0);
273: }
275: EXTERN_C_BEGIN
278: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscBool setmat)
279: {
280: ST_PRECOND *data = (ST_PRECOND*)st->data;
283: data->setmat = setmat;
284: return(0);
285: }
286: EXTERN_C_END
290: /*@
291: STPrecondSetKSPHasMat - Sets a flag indicating that during STSetUp the coefficient
292: matrix of the KSP linear system (A) must be set to be the same matrix as the
293: preconditioner (P).
295: Collective on ST
297: Input Parameter:
298: + st - the spectral transformation context
299: - setmat - the flag
301: Notes:
302: In most cases, the preconditioner matrix is used only in the PC object, but
303: in external solvers this matrix must be provided also as the A-matrix in
304: the KSP object.
306: Level: developer
308: .seealso: STPrecondGetKSPHasMat(), STSetShift()
309: @*/
310: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
311: {
317: PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
318: return(0);
319: }
321: EXTERN_C_BEGIN
324: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
325: {
326: ST_PRECOND *data = (ST_PRECOND*)st->data;
329: *setmat = data->setmat;
330: return(0);
331: }
332: EXTERN_C_END
336: /*@
337: STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
338: matrix of the KSP linear system (A) is set to be the same matrix as the
339: preconditioner (P).
341: Not Collective
343: Input Parameter:
344: . st - the spectral transformation context
346: Output Parameter:
347: . setmat - the flag
349: Level: developer
351: .seealso: STPrecondSetKSPHasMat(), STSetShift()
352: @*/
353: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
354: {
360: PetscTryMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
361: return(0);
362: }
366: PetscErrorCode STDestroy_Precond(ST st)
367: {
371: PetscFree(st->data);
372: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","",PETSC_NULL);
373: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","",PETSC_NULL);
374: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","",PETSC_NULL);
375: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","",PETSC_NULL);
376: return(0);
377: }
379: EXTERN_C_BEGIN
382: PetscErrorCode STCreate_Precond(ST st)
383: {
387: PetscNewLog(st,ST_PRECOND,&st->data);
388: st->ops->getbilinearform = STGetBilinearForm_Default;
389: st->ops->setup = STSetUp_Precond;
390: st->ops->setshift = STSetShift_Precond;
391: st->ops->destroy = STDestroy_Precond;
392: st->ops->setfromoptions = STSetFromOptions_Precond;
394: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","STPrecondGetMatForPC_Precond",STPrecondGetMatForPC_Precond);
395: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","STPrecondSetMatForPC_Precond",STPrecondSetMatForPC_Precond);
396: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","STPrecondGetKSPHasMat_Precond",STPrecondGetKSPHasMat_Precond);
397: PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","STPrecondSetKSPHasMat_Precond",STPrecondSetKSPHasMat_Precond);
399: STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
400: return(0);
401: }
402: EXTERN_C_END