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: VecScale(itg->xtilde[curl],1.0/norm);
121: itg->curl++;
122: }
123: return(0);
124: }
126: /* ---------------------------------------Method 2------------------------------------------------------------*/
127: typedef struct {
128: PetscInt method, /* 1 or 2 */
129: curl, /* Current number of basis vectors */
130: maxl, /* Maximum number of basis vectors */
131: refcnt;
132: PetscTruth monitor;
133: Mat mat;
134: KSP ksp;
135: PetscScalar *alpha; /* */
136: Vec *xtilde; /* Saved x vectors */
137: Vec Ax,guess;
138: } KSPFischerGuess_Method2;
142: PetscErrorCode KSPFischerGuessCreate_Method2(KSP ksp,int maxl,KSPFischerGuess_Method2 **ITG)
143: {
144: KSPFischerGuess_Method2 *itg;
145: PetscErrorCode ierr;
149: PetscMalloc(sizeof(KSPFischerGuess_Method2),&itg);
150: PetscMalloc(maxl * sizeof(PetscScalar),&itg->alpha);
151: PetscLogObjectMemory(ksp,sizeof(KSPFischerGuess_Method2) + maxl*sizeof(PetscScalar));
152: KSPGetVecs(ksp,maxl,&itg->xtilde,0,PETSC_NULL);
153: PetscLogObjectParents(ksp,maxl,itg->xtilde);
154: VecDuplicate(itg->xtilde[0],&itg->Ax);
155: PetscLogObjectParent(ksp,itg->Ax);
156: VecDuplicate(itg->xtilde[0],&itg->guess);
157: PetscLogObjectParent(ksp,itg->guess);
158: *ITG = itg;
159: return(0);
160: }
165: PetscErrorCode KSPFischerGuessDestroy_Method2(KSPFischerGuess_Method2 *itg)
166: {
170: PetscFree(itg->alpha);
171: VecDestroyVecs(itg->xtilde,itg->maxl);
172: VecDestroy(itg->Ax);
173: VecDestroy(itg->guess);
174: PetscFree(itg);
175: return(0);
176: }
179: /*
180: Given a basis generated already this computes a new guess x from the new right hand side b
181: Figures out the components of b in each btilde direction and adds them to x
182: */
185: PetscErrorCode KSPFischerGuessFormGuess_Method2(KSPFischerGuess_Method2 *itg,Vec b,Vec x)
186: {
188: PetscInt i;
193: VecSet(x,0.0);
194: VecMDot(b,itg->curl,itg->xtilde,itg->alpha);
195: if (itg->monitor) {
196: PetscPrintf(((PetscObject)itg->ksp)->comm,"KSPFischerGuess alphas = ");
197: for (i=0; i<itg->curl; i++ ){
198: PetscPrintf(((PetscObject)itg->ksp)->comm,"%G ",PetscAbsScalar(itg->alpha[i]));
199: }
200: PetscPrintf(((PetscObject)itg->ksp)->comm,"\n");
201: }
202: VecMAXPY(x,itg->curl,itg->alpha,itg->xtilde);
203: VecCopy(x,itg->guess);
204: /* Note: do not change the b right hand side as is done in the publication */
205: return(0);
206: }
210: PetscErrorCode KSPFischerGuessUpdate_Method2(KSPFischerGuess_Method2 *itg,Vec x)
211: {
212: PetscScalar norm;
214: int curl = itg->curl,i;
219: if (curl == itg->maxl) {
220: KSP_MatMult(itg->ksp,itg->mat,x,itg->Ax); /* norm = sqrt(x'Ax) */
221: VecDot(x,itg->Ax,&norm);
222: VecCopy(x,itg->xtilde[0]);
223: VecScale(itg->xtilde[0],1.0/PetscSqrtScalar(norm));
224: itg->curl = 1;
225: } else {
226: if (!curl) {
227: VecCopy(x,itg->xtilde[curl]);
228: } else {
229: VecWAXPY(itg->xtilde[curl],-1.0,itg->guess,x);
230: }
231: KSP_MatMult(itg->ksp,itg->mat,itg->xtilde[curl],itg->Ax);
232: VecMDot(itg->Ax,curl,itg->xtilde,itg->alpha);
233: for (i=0; i<curl; i++) itg->alpha[i] = -itg->alpha[i];
234: VecMAXPY(itg->xtilde[curl],curl,itg->alpha,itg->xtilde);
236: KSP_MatMult(itg->ksp,itg->mat,itg->xtilde[curl],itg->Ax); /* norm = sqrt(xtilde[curl]'Axtilde[curl]) */
237: VecDot(itg->xtilde[curl],itg->Ax,&norm);
238: VecScale(itg->xtilde[curl],1.0/PetscSqrtScalar(norm));
239: itg->curl++;
240: }
241: return(0);
242: }
244: /* ---------------------------------------------------------------------------------------------------------*/
247: /*@C
248: KSPFischerGuessCreate - Implements Paul Fischer's initial guess algorithm Method 1 and 2 for situations where
249: a linear system is solved repeatedly
251: References:
252: http://ntrs.nasa.gov/archive/nasa/casi.ntrs.nasa.gov/19940020363_1994020363.pdf
254: Notes: the algorithm is different from the paper because we do not CHANGE the right hand side of the new
255: problem and solve the problem with an initial guess of zero, rather we solve the original new problem
256: with a nonzero initial guess (this is done so that the linear solver convergence tests are based on
257: the original RHS.) But we use the xtilde = x - xguess as the new direction so that it is not
258: mostly orthogonal to the previous solutions.
260: These are not intended to be used directly, they are called by KSP automatically when the
261: KSP option KSPSetFischerGuess(KSP,PetscInt,PetscInt) or -ksp_guess_fischer <int,int>
263: Method 2 is only for positive definite matrices, since it uses the A norm.
265: This is not currently programmed as a PETSc class because there are only two methods; if more methods
266: are introduced it should be changed. For example the Knoll guess should be included
268: Level: advanced
270: @*/
271: PetscErrorCode KSPFischerGuessCreate(KSP ksp,PetscInt method,PetscInt maxl,KSPFischerGuess *itg)
272: {
273: PetscErrorCode ierr;
276: if (method == 1) {
277: KSPFischerGuessCreate_Method1(ksp,maxl,(KSPFischerGuess_Method1 **)itg);
278: } else if (method == 2) {
279: KSPFischerGuessCreate_Method2(ksp,maxl,(KSPFischerGuess_Method2 **)itg);
280: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
281: (*itg)->method = method;
282: (*itg)->curl = 0;
283: (*itg)->maxl = maxl;
284: (*itg)->ksp = ksp;
285: (*itg)->refcnt = 1;
286: KSPGetOperators(ksp,&(*itg)->mat,PETSC_NULL,PETSC_NULL);
287: return(0);
288: }
292: PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess ITG)
293: {
294: PetscErrorCode ierr;
297: PetscOptionsGetTruth(((PetscObject)ITG->ksp)->prefix,"-ksp_fischer_guess_monitor",&ITG->monitor,PETSC_NULL);
298: return(0);
299: }
303: PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess ITG)
304: {
305: PetscErrorCode ierr;
308: if (--ITG->refcnt) return(0);
310: if (ITG->method == 1) {
311: KSPFischerGuessDestroy_Method1((KSPFischerGuess_Method1 *)ITG);
312: } else if (ITG->method == 2) {
313: KSPFischerGuessDestroy_Method2((KSPFischerGuess_Method2 *)ITG);
314: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
315: return(0);
316: }
320: PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess itg,Vec x)
321: {
322: PetscErrorCode ierr;
325: if (itg->method == 1) {
326: KSPFischerGuessUpdate_Method1((KSPFischerGuess_Method1 *)itg,x);
327: } else if (itg->method == 2) {
328: KSPFischerGuessUpdate_Method2((KSPFischerGuess_Method2 *)itg,x);
329: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
330: return(0);
331: }
335: PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess itg,Vec b,Vec x)
336: {
337: PetscErrorCode ierr;
340: if (itg->method == 1) {
341: KSPFischerGuessFormGuess_Method1((KSPFischerGuess_Method1 *)itg,b,x);
342: } else if (itg->method == 2) {
343: KSPFischerGuessFormGuess_Method2((KSPFischerGuess_Method2 *)itg,b,x);
344: } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Method can only be 1 or 2");
345: return(0);
346: }
350: /*
351: KSPFischerGuessReset - This is called whenever KSPSetOperators() is called to tell the
352: initial guess object that the matrix is changed and so the initial guess object
353: must restart from scratch building the subspace where the guess is computed from.
354: */
355: PetscErrorCode KSPFischerGuessReset(KSPFischerGuess itg)
356: {
360: itg->curl = 0;
361: KSPGetOperators(itg->ksp,&itg->mat,PETSC_NULL,PETSC_NULL);
362: return(0);
363: }