Actual source code: iguess.c

 2:  #include src/ksp/ksp/kspimpl.h
  3: /* 
  4:   This code inplements Paul Fischer's initial guess code for situations where
  5:   a linear system is solved repeatedly 
  6:  */

  8: typedef struct {
  9:     int         curl,     /* Current number of basis vectors */
 10:                 maxl;     /* Maximum number of basis vectors */
 11:     PetscScalar *alpha;   /* */
 12:     Vec         *xtilde,  /* Saved x vectors */
 13:                 *btilde;  /* Saved b vectors */
 14: } KSPIGUESS;

 18: PetscErrorCode KSPGuessCreate(KSP ksp,int  maxl,void **ITG)
 19: {
 20:   KSPIGUESS *itg;

 23:   *ITG = 0;
 26:   PetscMalloc(sizeof(KSPIGUESS),&itg);
 27:   itg->curl = 0;
 28:   itg->maxl = maxl;
 29:   PetscMalloc(maxl * sizeof(PetscScalar),&itg->alpha);
 30:   PetscLogObjectMemory(ksp,sizeof(KSPIGUESS) + maxl*sizeof(PetscScalar));
 31:   KSPGetVecs(ksp,maxl,&itg->xtilde);
 32:   PetscLogObjectParents(ksp,maxl,itg->xtilde);
 33:   KSPGetVecs(ksp,maxl,&itg->btilde);
 34:   PetscLogObjectParents(ksp,maxl,itg->btilde);
 35:   *ITG = (void*)itg;
 36:   return(0);
 37: }

 41: PetscErrorCode KSPGuessDestroy(KSP ksp,KSPIGUESS *itg)
 42: {

 47:   PetscFree(itg->alpha);
 48:   VecDestroyVecs(itg->btilde,itg->maxl);
 49:   VecDestroyVecs(itg->xtilde,itg->maxl);
 50:   PetscFree(itg);
 51:   return(0);
 52: }

 56: PetscErrorCode KSPGuessFormB(KSP ksp,KSPIGUESS *itg,Vec b)
 57: {
 59:   int         i;
 60:   PetscScalar tmp;

 66:   for (i=1; i<=itg->curl; i++) {
 67:     VecDot(itg->btilde[i-1],b,&(itg->alpha[i-1]));
 68:     tmp = -itg->alpha[i-1];
 69:     VecAXPY(&tmp,itg->btilde[i-1],b);
 70:   }
 71:   return(0);
 72: }

 76: PetscErrorCode KSPGuessFormX(KSP ksp,KSPIGUESS *itg,Vec x)
 77: {
 79:   int i;

 85:   VecCopy(x,itg->xtilde[itg->curl]);
 86:   for (i=1; i<=itg->curl; i++) {
 87:     VecAXPY(&itg->alpha[i-1],itg->xtilde[i-1],x);
 88:   }
 89:   return(0);
 90: }

 94: PetscErrorCode  KSPGuessUpdate(KSP ksp,Vec x,KSPIGUESS *itg)
 95: {
 96:   PetscReal    normax,norm;
 97:   PetscScalar  tmp;
 98:   MatStructure pflag;
100:   int          curl = itg->curl,i;
101:   Mat          Amat,Pmat;

107:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
108:   if (curl == itg->maxl) {
109:     KSP_MatMult(ksp,Amat,x,itg->btilde[0]);
110:     VecNorm(itg->btilde[0],NORM_2,&normax);
111:     tmp = 1.0/normax; VecScale(&tmp,itg->btilde[0]);
112:     /* VCOPY(ksp->vc,x,itg->xtilde[0]); */
113:     VecScale(&tmp,itg->xtilde[0]);
114:   } else {
115:     KSP_MatMult(ksp,Amat,itg->xtilde[curl],itg->btilde[curl]);
116:     for (i=1; i<=curl; i++) {
117:       VecDot(itg->btilde[curl],itg->btilde[i-1],itg->alpha+i-1);
118:     }
119:     for (i=1; i<=curl; i++) {
120:       tmp  = -itg->alpha[i-1];
121:       VecAXPY(&tmp,itg->btilde[i-1],itg->btilde[curl]);
122:       VecAXPY(&itg->alpha[i-1],itg->xtilde[i-1],itg->xtilde[curl]);
123:     }
124:     VecNorm(itg->btilde[curl],NORM_2,&norm);
125:     tmp = 1.0/norm; VecScale(&tmp,itg->btilde[curl]);
126:     VecNorm(itg->xtilde[curl],NORM_2,&norm);
127:     VecScale(&tmp,itg->xtilde[curl]);
128:     itg->curl++;
129:   }
130:   return(0);
131: }