Actual source code: iguess.c

  1: /*$Id: iguess.c,v 1.38 2001/08/07 03:03:45 balay Exp $*/

 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: int KSPGuessCreate(KSP ksp,int  maxl,void **ITG)
 20: {
 21:   KSPIGUESS *itg;
 22:   int       ierr;

 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:   VecDuplicateVecs(ksp->vec_rhs,maxl,&itg->xtilde);
 33:   PetscLogObjectParents(ksp,maxl,itg->xtilde);
 34:   VecDuplicateVecs(ksp->vec_rhs,maxl,&itg->btilde);
 35:   PetscLogObjectParents(ksp,maxl,itg->btilde);
 36:   *ITG = (void *)itg;
 37:   return(0);
 38: }

 42: int 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: int KSPGuessFormB(KSP ksp,KSPIGUESS *itg,Vec b)
 58: {
 59:   int         i,ierr;
 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: int KSPGuessFormX(KSP ksp,KSPIGUESS *itg,Vec x)
 77: {
 78:   int i,ierr;

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

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

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