Actual source code: iguess.c
1: #define PETSCKSP_DLL
3: #include src/ksp/ksp/kspimpl.h
4: /*
5: This code inplements Paul Fischer's initial guess code for situations where
6: a linear system is solved repeatedly
7: */
9: typedef struct {
10: int curl, /* Current number of basis vectors */
11: maxl; /* Maximum number of basis vectors */
12: PetscScalar *alpha; /* */
13: Vec *xtilde, /* Saved x vectors */
14: *btilde; /* Saved b vectors */
15: } KSPIGUESS;
19: PetscErrorCode PETSCKSP_DLLEXPORT KSPGuessCreate(KSP ksp,int maxl,void **ITG)
20: {
21: KSPIGUESS *itg;
24: *ITG = 0;
27: PetscMalloc(sizeof(KSPIGUESS),&itg);
28: itg->curl = 0;
29: itg->maxl = maxl;
30: PetscMalloc(maxl * sizeof(PetscScalar),&itg->alpha);
31: PetscLogObjectMemory(ksp,sizeof(KSPIGUESS) + maxl*sizeof(PetscScalar));
32: KSPGetVecs(ksp,maxl,&itg->xtilde,0,PETSC_NULL);
33: PetscLogObjectParents(ksp,maxl,itg->xtilde);
34: KSPGetVecs(ksp,maxl,&itg->btilde,0,PETSC_NULL);
35: PetscLogObjectParents(ksp,maxl,itg->btilde);
36: *ITG = (void*)itg;
37: return(0);
38: }
42: PetscErrorCode PETSCKSP_DLLEXPORT KSPGuessDestroy(KSP ksp,KSPIGUESS *itg)
43: {
48: PetscFree(itg->alpha);
49: VecDestroyVecs(itg->btilde,itg->maxl);
50: VecDestroyVecs(itg->xtilde,itg->maxl);
51: PetscFree(itg);
52: return(0);
53: }
57: PetscErrorCode PETSCKSP_DLLEXPORT KSPGuessFormB(KSP ksp,KSPIGUESS *itg,Vec b)
58: {
60: PetscInt i;
66: for (i=1; i<=itg->curl; i++) {
67: VecDot(itg->btilde[i-1],b,&(itg->alpha[i-1]));
68: VecAXPY(b,-itg->alpha[i-1],itg->btilde[i-1]);
69: }
70: return(0);
71: }
75: PetscErrorCode PETSCKSP_DLLEXPORT KSPGuessFormX(KSP ksp,KSPIGUESS *itg,Vec x)
76: {
78: int i;
84: VecCopy(x,itg->xtilde[itg->curl]);
85: for (i=1; i<=itg->curl; i++) {
86: VecAXPY(x,itg->alpha[i-1],itg->xtilde[i-1]);
87: }
88: return(0);
89: }
93: PetscErrorCode PETSCKSP_DLLEXPORT KSPGuessUpdate(KSP ksp,Vec x,KSPIGUESS *itg)
94: {
95: PetscReal normax,norm;
96: PetscScalar tmp;
97: MatStructure pflag;
99: int curl = itg->curl,i;
100: Mat Amat,Pmat;
106: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
107: if (curl == itg->maxl) {
108: KSP_MatMult(ksp,Amat,x,itg->btilde[0]);
109: VecNorm(itg->btilde[0],NORM_2,&normax);
110: tmp = 1.0/normax; VecScale(itg->btilde[0],tmp);
111: /* VCOPY(ksp->vc,x,itg->xtilde[0]); */
112: VecScale(itg->xtilde[0],tmp);
113: } else {
114: KSP_MatMult(ksp,Amat,itg->xtilde[curl],itg->btilde[curl]);
115: for (i=1; i<=curl; i++) {
116: VecDot(itg->btilde[curl],itg->btilde[i-1],itg->alpha+i-1);
117: }
118: for (i=1; i<=curl; i++) {
119: VecAXPY(itg->btilde[curl],-itg->alpha[i-1],itg->btilde[i-1]);
120: VecAXPY(itg->xtilde[curl],itg->alpha[i-1],itg->xtilde[i-1]);
121: }
122: VecNormalize(itg->btilde[curl],&norm);
123: VecNormalize(itg->xtilde[curl],&norm);
124: itg->curl++;
125: }
126: return(0);
127: }