Actual source code: tfqmr.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/kspimpl.h>

  6: static PetscErrorCode KSPSetUp_TFQMR(KSP ksp)
  7: {

 11:   if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"no symmetric preconditioning for KSPTFQMR");
 12:   KSPSetWorkVecs(ksp,9);
 13:   return(0);
 14: }

 18: static PetscErrorCode  KSPSolve_TFQMR(KSP ksp)
 19: {
 21:   PetscInt       i,m;
 22:   PetscScalar    rho,rhoold,a,s,b,eta,etaold,psiold,cf;
 23:   PetscReal      dp,dpold,w,dpest,tau,psi,cm;
 24:   Vec            X,B,V,P,R,RP,T,T1,Q,U,D,AUQ;

 27:   X   = ksp->vec_sol;
 28:   B   = ksp->vec_rhs;
 29:   R   = ksp->work[0];
 30:   RP  = ksp->work[1];
 31:   V   = ksp->work[2];
 32:   T   = ksp->work[3];
 33:   Q   = ksp->work[4];
 34:   P   = ksp->work[5];
 35:   U   = ksp->work[6];
 36:   D   = ksp->work[7];
 37:   T1  = ksp->work[8];
 38:   AUQ = V;

 40:   /* Compute initial preconditioned residual */
 41:   KSPInitialResidual(ksp,X,V,T,R,B);

 43:   /* Test for nothing to do */
 44:   VecNorm(R,NORM_2,&dp);
 45:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 46:   ksp->rnorm = dp;
 47:   ksp->its   = 0;
 48:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 49:   KSPMonitor(ksp,0,dp);
 50:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 51:   if (ksp->reason) return(0);

 53:   /* Make the initial Rp == R */
 54:   VecCopy(R,RP);

 56:   /* Set the initial conditions */
 57:   etaold = 0.0;
 58:   psiold = 0.0;
 59:   tau    = dp;
 60:   dpold  = dp;

 62:   VecDot(R,RP,&rhoold);       /* rhoold = (r,rp)     */
 63:   VecCopy(R,U);
 64:   VecCopy(R,P);
 65:   KSP_PCApplyBAorAB(ksp,P,V,T);
 66:   VecSet(D,0.0);

 68:   i=0;
 69:   do {
 70:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
 71:     ksp->its++;
 72:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
 73:     VecDot(V,RP,&s);          /* s <- (v,rp)          */
 74:     a    = rhoold / s;                              /* a <- rho / s         */
 75:     VecWAXPY(Q,-a,V,U);  /* q <- u - a v         */
 76:     VecWAXPY(T,1.0,U,Q);     /* t <- u + q           */
 77:     KSP_PCApplyBAorAB(ksp,T,AUQ,T1);
 78:     VecAXPY(R,-a,AUQ);      /* r <- r - a K (u + q) */
 79:     VecNorm(R,NORM_2,&dp);
 80:     for (m=0; m<2; m++) {
 81:       if (!m) w = PetscSqrtReal(dp*dpold);
 82:       else w = dp;
 83:       psi = w / tau;
 84:       cm  = 1.0 / PetscSqrtReal(1.0 + psi * psi);
 85:       tau = tau * psi * cm;
 86:       eta = cm * cm * a;
 87:       cf  = psiold * psiold * etaold / a;
 88:       if (!m) {
 89:         VecAYPX(D,cf,U);
 90:       } else {
 91:         VecAYPX(D,cf,Q);
 92:       }
 93:       VecAXPY(X,eta,D);

 95:       dpest      = PetscSqrtReal(m + 1.0) * tau;
 96:       PetscObjectSAWsTakeAccess((PetscObject)ksp);
 97:       ksp->rnorm = dpest;
 98:       PetscObjectSAWsGrantAccess((PetscObject)ksp);
 99:       KSPLogResidualHistory(ksp,dpest);
100:       KSPMonitor(ksp,i+1,dpest);
101:       (*ksp->converged)(ksp,i+1,dpest,&ksp->reason,ksp->cnvP);
102:       if (ksp->reason) break;

104:       etaold = eta;
105:       psiold = psi;
106:     }
107:     if (ksp->reason) break;

109:     VecDot(R,RP,&rho);        /* rho <- (r,rp)       */
110:     b    = rho / rhoold;                            /* b <- rho / rhoold   */
111:     VecWAXPY(U,b,Q,R);       /* u <- r + b q        */
112:     VecAXPY(Q,b,P);
113:     VecWAXPY(P,b,Q,U);       /* p <- u + b(q + b p) */
114:     KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p  */

116:     rhoold = rho;
117:     dpold  = dp;

119:     i++;
120:   } while (i<ksp->max_it);
121:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

123:   KSPUnwindPreconditioner(ksp,X,T);
124:   return(0);
125: }

127: /*MC
128:      KSPTFQMR - A transpose free QMR (quasi minimal residual),

130:    Options Database Keys:
131: .   see KSPSolve()

133:    Level: beginner

135:    Notes: Supports left and right preconditioning, but not symmetric

137:           The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
138:           That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
139:           it is a bound on the true residual.

141:    References: Freund, 1993

143: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTCQMR
144: M*/
147: PETSC_EXTERN PetscErrorCode KSPCreate_TFQMR(KSP ksp)
148: {

152:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
153:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);

155:   ksp->data                = (void*)0;
156:   ksp->ops->setup          = KSPSetUp_TFQMR;
157:   ksp->ops->solve          = KSPSolve_TFQMR;
158:   ksp->ops->destroy        = KSPDestroyDefault;
159:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
160:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
161:   ksp->ops->setfromoptions = 0;
162:   ksp->ops->view           = 0;
163:   return(0);
164: }