Actual source code: iguess.c

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

 3:  #include src/sles/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;

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

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

 38: int KSPGuessDestroy(KSP ksp,KSPIGUESS *itg)
 39: {

 44:   PetscFree(itg->alpha);
 45:   VecDestroyVecs(itg->btilde,itg->maxl);
 46:   VecDestroyVecs(itg->xtilde,itg->maxl);
 47:   PetscFree(itg);
 48:   return(0);
 49: }

 51: int KSPGuessFormB(KSP ksp,KSPIGUESS *itg,Vec b)
 52: {
 53:   int    i,ierr;
 54:   PetscScalar tmp;

 58:   for (i=1; i<=itg->curl; i++) {
 59:     VecDot(itg->btilde[i-1],b,&(itg->alpha[i-1]));
 60:     tmp = -itg->alpha[i-1];
 61:     VecAXPY(&tmp,itg->btilde[i-1],b);
 62:   }
 63:   return(0);
 64: }

 66: int KSPGuessFormX(KSP ksp,KSPIGUESS *itg,Vec x)
 67: {
 68:   int i,ierr;

 72:   VecCopy(x,itg->xtilde[itg->curl]);
 73:   for (i=1; i<=itg->curl; i++) {
 74:     VecAXPY(&itg->alpha[i-1],itg->xtilde[i-1],x);
 75:   }
 76:   return(0);
 77: }

 79: int  KSPGuessUpdate(KSP ksp,Vec x,KSPIGUESS *itg)
 80: {
 81:   PetscReal    normax,norm;
 82:   PetscScalar  tmp;
 83:   MatStructure pflag;
 84:   int          curl = itg->curl,i,ierr;
 85:   Mat          Amat,Pmat;

 89:   PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
 90:   if (curl == itg->maxl) {
 91:     KSP_MatMult(ksp,Amat,x,itg->btilde[0]);
 92:     VecNorm(itg->btilde[0],NORM_2,&normax);
 93:     tmp = 1.0/normax; VecScale(&tmp,itg->btilde[0]);
 94:     /* VCOPY(ksp->vc,x,itg->xtilde[0]); */
 95:     VecScale(&tmp,itg->xtilde[0]);
 96:   } else {
 97:     KSP_MatMult(ksp,Amat,itg->xtilde[curl],itg->btilde[curl]);
 98:     for (i=1; i<=curl; i++) {
 99:       VecDot(itg->btilde[curl],itg->btilde[i-1],itg->alpha+i-1);
100:     }
101:     for (i=1; i<=curl; i++) {
102:       tmp  = -itg->alpha[i-1];
103:       VecAXPY(&tmp,itg->btilde[i-1],itg->btilde[curl]);
104:       VecAXPY(&itg->alpha[i-1],itg->xtilde[i-1],itg->xtilde[curl]);
105:     }
106:     VecNorm(itg->btilde[curl],NORM_2,&norm);
107:     tmp = 1.0/norm; VecScale(&tmp,itg->btilde[curl]);
108:     VecNorm(itg->xtilde[curl],NORM_2,&norm);
109:     VecScale(&tmp,itg->xtilde[curl]);
110:     itg->curl++;
111:   }
112:   return(0);
113: }