Actual source code: iguess.c
1: #define PETSCKSP_DLL
3: #include private/kspimpl.h
5: /* ---------------------------------------Method 1------------------------------------------------------------*/
6: typedef struct {
7: PetscInt method, /* 1 or 2 */
8: curl, /* Current number of basis vectors */
9: maxl, /* Maximum number of basis vectors */
10: refcnt;
11: PetscTruth monitor;
12: Mat mat;
13: KSP ksp;
14: PetscScalar *alpha; /* */
15: Vec *xtilde, /* Saved x vectors */
16: *btilde; /* Saved b vectors */
17: Vec guess;
18: } KSPFischerGuess_Method1;
23: PetscErrorCode KSPFischerGuessCreate_Method1(KSP ksp,int maxl,KSPFischerGuess_Method1 **ITG)
24: {
25: KSPFischerGuess_Method1 *itg;
26: PetscErrorCode ierr;
30: PetscMalloc(sizeof(KSPFischerGuess_Method1),&itg);
31: PetscMalloc(maxl * sizeof(PetscScalar),&itg->alpha);
32: PetscLogObjectMemory(ksp,sizeof(KSPFischerGuess_Method1) + maxl*sizeof(PetscScalar));
33: KSPGetVecs(ksp,maxl,&itg->xtilde,0,PETSC_NULL);
34: PetscLogObjectParents(ksp,maxl,itg->xtilde);
35: KSPGetVecs(ksp,maxl,&itg->btilde,0,PETSC_NULL);
36: PetscLogObjectParents(ksp,maxl,itg->btilde);
37: VecDuplicate(itg->xtilde[0],&itg->guess);
38: PetscLogObjectParent(ksp,itg->guess);
39: *ITG = itg;
40: return(0);
41: }
46: PetscErrorCode KSPFischerGuessDestroy_Method1(KSPFischerGuess_Method1 *itg)
47: {
51: PetscFree(itg->alpha);
52: VecDestroyVecs(itg->btilde,itg->maxl);
53: VecDestroyVecs(itg->xtilde,itg->maxl);
54: VecDestroy(itg->guess);
55: PetscFree(itg);
56: return(0);
57: }
60: /*
61: Given a basis generated already this computes a new guess x from the new right hand side b
62: Figures out the components of b in each btilde direction and adds them to x
63: */
66: PetscErrorCode KSPFischerGuessFormGuess_Method1(KSPFischerGuess_Method1 *itg,Vec b,Vec x)
67: {
69: PetscInt i;
74: VecSet(x,0.0);
75: VecMDot(b,itg->curl,itg->btilde,itg->alpha);
76: if (itg->monitor) {
77: PetscPrintf(((PetscObject)itg->ksp)->comm,"KSPFischerGuess alphas = ");
78: for (i=0; i<itg->curl; i++ ){
79: PetscPrintf(((PetscObject)itg->ksp)->comm,"%G ",PetscAbsScalar(itg->alpha[i]));
80: }
81: PetscPrintf(((PetscObject)itg->ksp)->comm,"\n");
82: }
83: VecMAXPY(x,itg->curl,itg->alpha,itg->xtilde);
84: VecCopy(x,itg->guess);
85: /* Note: do not change the b right hand side as is done in the publication */
86: return(0);
87: }
91: PetscErrorCode KSPFischerGuessUpdate_Method1(KSPFischerGuess_Method1 *itg,Vec x)
92: {
93: PetscReal norm;
95: int curl = itg->curl,i;
100: if (curl == itg->maxl) {
101: KSP_MatMult(itg->ksp,itg->mat,x,itg->btilde[0]);
102: VecNormalize(itg->btilde[0],&norm);
103: VecCopy(x,itg->xtilde[0]);
104: VecScale(itg->xtilde[0],1.0/norm);
105: itg->curl = 1;
106: } else {
107: if (!curl) {
108: VecCopy(x,itg->xtilde[curl]);
109: } else {
110: VecWAXPY(itg->xtilde[curl],-1.0,itg->guess,x);
111: }
113: KSP_MatMult(itg->ksp,itg->mat,itg->xtilde[curl],itg->btilde[curl]);
114: VecMDot(itg->btilde[curl],curl,itg->btilde,itg->alpha);
115: for (i=0; i<curl; i++) itg->alpha[i] = -itg->alpha[i];
116: VecMAXPY(itg->btilde[curl],curl,itg->alpha,itg->btilde);
117: VecMAXPY(itg->xtilde[curl],curl,itg->alpha,itg->xtilde);
119: VecNormalize(itg->btilde[curl],&norm);
120: if (norm) {
121: VecScale(itg->xtilde[curl],1.0/norm);
122: itg->curl++;
123: } else {
124: PetscInfo(itg->ksp,"Not increasing dimension of Fischer space because new direction is identical to previous");
125: }
126: }
127: return(0);
128: }
130: /* ---------------------------------------Method 2------------------------------------------------------------*/
131: typedef struct {
132: PetscInt method, /* 1 or 2 */
133: curl, /* Current number of basis vectors */
134: maxl, /* Maximum number of basis vectors */
135: refcnt;
136: PetscTruth monitor;
137: Mat mat;
138: KSP ksp;
139: PetscScalar *alpha; /* */
140: Vec *xtilde; /* Saved x vectors */
141: Vec Ax,guess;
142: } KSPFischerGuess_Method2;
146: PetscErrorCode KSPFischerGuessCreate_Method2(KSP ksp,int maxl,KSPFischerGuess_Method2 **ITG)
147: {
148: KSPFischerGuess_Method2 *itg;
149: PetscErrorCode ierr;
153: PetscMalloc(sizeof(KSPFischerGuess_Method2),&itg);
154: PetscMalloc(maxl * sizeof(PetscScalar),&itg->alpha);
155: PetscLogObjectMemory(ksp,sizeof(KSPFischerGuess_Method2) + maxl*sizeof(PetscScalar));
156: KSPGetVecs(ksp,maxl,&itg->xtilde,0,PETSC_NULL);
157: PetscLogObjectParents(ksp,maxl,itg->xtilde);
158: VecDuplicate(itg->xtilde[0],&itg->Ax);
159: PetscLogObjectParent(ksp,itg->Ax);
160: VecDuplicate(itg->xtilde[0],&itg->guess);
161: PetscLogObjectParent(ksp,itg->guess);
162: *ITG = itg;
163: return(0);
164: }
169: PetscErrorCode KSPFischerGuessDestroy_Method2(KSPFischerGuess_Method2 *itg)
170: {
174: PetscFree(itg->alpha);
175: VecDestroyVecs(itg->xtilde,itg->maxl);
176: VecDestroy(itg->Ax);
177: VecDestroy(itg->guess);
178: PetscFree(itg);
179: return(0);
180: }
183: /*
184: Given a basis generated already this computes a new guess x from the new right hand side b
185: Figures out the components of b in each btilde direction and adds them to x
186: */
189: PetscErrorCode KSPFischerGuessFormGuess_Method2(KSPFischerGuess_Method2 *itg,Vec b,Vec x)
190: {
192: PetscInt i;
197: VecSet(x,0.0);
198: VecMDot(b,itg->curl,itg->xtilde,itg->alpha);
199: if (itg->monitor) {
200: PetscPrintf(((PetscObject)itg->ksp)->comm,"KSPFischerGuess alphas = ");
201: for (i=0; i<itg->curl; i++ ){
202: PetscPrintf(((PetscObject)itg->ksp)->comm,"%G ",PetscAbsScalar(itg->alpha[i]));
203: }
204: PetscPrintf(((PetscObject)itg->ksp)->comm,"\n");
205: }
206: VecMAXPY(x,itg->curl,itg->alpha,itg->xtilde);
207: VecCopy(x,itg->guess);
208: /* Note: do not change the b right hand side as is done in the publication */
209: return(0);
210: }
214: PetscErrorCode KSPFischerGuessUpdate_Method2(KSPFischerGuess_Method2 *itg,Vec x)
215: {
216: PetscScalar norm;
218: int curl = itg->curl,i;
223: if (curl == itg->maxl) {
224: KSP_MatMult(itg->ksp,itg->mat,x,itg->Ax); /* norm = sqrt(x'Ax) */
225: VecDot(x,itg->Ax,&norm);
226: VecCopy(x,itg->xtilde[0]);
227: VecScale(itg->xtilde[0],1.0/PetscSqrtScalar(norm));
228: itg->curl = 1;
229: } else {
230: if (!curl) {
231: VecCopy(x,itg->xtilde[curl]);
232: } else {
233: VecWAXPY(itg->xtilde[curl],-1.0,itg->guess,x);
234: }
235: KSP_MatMult(itg->ksp,itg->mat,itg->xtilde[curl],itg->Ax);
236: VecMDot(itg->Ax,curl,itg->xtilde,itg->alpha);
237: for (i=0; i<curl; i++) itg->alpha[i] = -itg->alpha[i];
238: VecMAXPY(itg->xtilde[curl],curl,itg->alpha,itg->xtilde);
240: KSP_MatMult(itg->ksp,itg->mat,itg->xtilde[curl],itg->Ax); /* norm = sqrt(xtilde[curl]'Axtilde[curl]) */
241: VecDot(itg->xtilde[curl],itg->Ax,&norm);
242: if (PetscAbsScalar(norm) != 0.0) {
243: VecScale(itg->xtilde[curl],1.0/PetscSqrtScalar(norm));
244: itg->curl++;
245: } else {
246: PetscInfo(itg->ksp,"Not increasing dimension of Fischer space because new direction is identical to previous");
247: }
248: }
249: return(0);
250: }
252: /* ---------------------------------------------------------------------------------------------------------*/
255: /*@C
256: KSPFischerGuessCreate - Implements Paul Fischer's initial guess algorithm Method 1 and 2 for situations where
257: a linear system is solved repeatedly
259: References:
260: http://ntrs.nasa.gov/archive/nasa/casi.ntrs.nasa.gov/19940020363_1994020363.pdf
262: Notes: the algorithm is different from the paper because we do not CHANGE the right hand side of the new
263: problem and solve the problem with an initial guess of zero, rather we solve the original new problem
264: with a nonzero initial guess (this is done so that the linear solver convergence tests are based on
265: the original RHS.) But we use the xtilde = x - xguess as the new direction so that it is not
266: mostly orthogonal to the previous solutions.
268: These are not intended to be used directly, they are called by KSP automatically when the
269: KSP option KSPSetFischerGuess(KSP,PetscInt,PetscInt) or -ksp_guess_fischer <int,int>
271: Method 2 is only for positive definite matrices, since it uses the A norm.
273: This is not currently programmed as a PETSc class because there are only two methods; if more methods
274: are introduced it should be changed. For example the Knoll guess should be included
276: Level: advanced
278: @*/
279: PetscErrorCode KSPFischerGuessCreate(KSP ksp,PetscInt method,PetscInt maxl,KSPFischerGuess *itg)
280: {
281: PetscErrorCode ierr;
284: if (method == 1) {
285: KSPFischerGuessCreate_Method1(ksp,maxl,(KSPFischerGuess_Method1 **)itg);
286: } else if (method == 2) {
287: KSPFischerGuessCreate_Method2(ksp,maxl,(KSPFischerGuess_Method2 **)itg);
288: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
289: (*itg)->method = method;
290: (*itg)->curl = 0;
291: (*itg)->maxl = maxl;
292: (*itg)->ksp = ksp;
293: (*itg)->refcnt = 1;
294: KSPGetOperators(ksp,&(*itg)->mat,PETSC_NULL,PETSC_NULL);
295: return(0);
296: }
300: PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess ITG)
301: {
302: PetscErrorCode ierr;
305: PetscOptionsGetTruth(((PetscObject)ITG->ksp)->prefix,"-ksp_fischer_guess_monitor",&ITG->monitor,PETSC_NULL);
306: return(0);
307: }
311: PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess ITG)
312: {
313: PetscErrorCode ierr;
316: if (--ITG->refcnt) return(0);
318: if (ITG->method == 1) {
319: KSPFischerGuessDestroy_Method1((KSPFischerGuess_Method1 *)ITG);
320: } else if (ITG->method == 2) {
321: KSPFischerGuessDestroy_Method2((KSPFischerGuess_Method2 *)ITG);
322: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
323: return(0);
324: }
328: PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess itg,Vec x)
329: {
330: PetscErrorCode ierr;
333: if (itg->method == 1) {
334: KSPFischerGuessUpdate_Method1((KSPFischerGuess_Method1 *)itg,x);
335: } else if (itg->method == 2) {
336: KSPFischerGuessUpdate_Method2((KSPFischerGuess_Method2 *)itg,x);
337: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
338: return(0);
339: }
343: PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess itg,Vec b,Vec x)
344: {
345: PetscErrorCode ierr;
348: if (itg->method == 1) {
349: KSPFischerGuessFormGuess_Method1((KSPFischerGuess_Method1 *)itg,b,x);
350: } else if (itg->method == 2) {
351: KSPFischerGuessFormGuess_Method2((KSPFischerGuess_Method2 *)itg,b,x);
352: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
353: return(0);
354: }
358: /*
359: KSPFischerGuessReset - This is called whenever KSPSetOperators() is called to tell the
360: initial guess object that the matrix is changed and so the initial guess object
361: must restart from scratch building the subspace where the guess is computed from.
362: */
363: PetscErrorCode KSPFischerGuessReset(KSPFischerGuess itg)
364: {
368: itg->curl = 0;
369: KSPGetOperators(itg->ksp,&itg->mat,PETSC_NULL,PETSC_NULL);
370: return(0);
371: }