Actual source code: tfqmr.c


  2: #include <petsc/private/kspimpl.h>

  4: static PetscErrorCode KSPSetUp_TFQMR(KSP ksp)
  5: {
  7:   KSPSetWorkVecs(ksp, 9);
  8:   return 0;
  9: }

 11: static PetscErrorCode KSPSolve_TFQMR(KSP ksp)
 12: {
 13:   PetscInt    i, m;
 14:   PetscScalar rho, rhoold, a, s, b, eta, etaold, psiold, cf;
 15:   PetscReal   dp, dpold, w, dpest, tau, psi, cm;
 16:   Vec         X, B, V, P, R, RP, T, T1, Q, U, D, AUQ;

 18:   X   = ksp->vec_sol;
 19:   B   = ksp->vec_rhs;
 20:   R   = ksp->work[0];
 21:   RP  = ksp->work[1];
 22:   V   = ksp->work[2];
 23:   T   = ksp->work[3];
 24:   Q   = ksp->work[4];
 25:   P   = ksp->work[5];
 26:   U   = ksp->work[6];
 27:   D   = ksp->work[7];
 28:   T1  = ksp->work[8];
 29:   AUQ = V;

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

 34:   /* Test for nothing to do */
 35:   VecNorm(R, NORM_2, &dp);
 36:   KSPCheckNorm(ksp, dp);
 37:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 38:   if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dp;
 39:   else ksp->rnorm = 0.0;
 40:   ksp->its = 0;
 41:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 42:   KSPMonitor(ksp, 0, ksp->rnorm);
 43:   (*ksp->converged)(ksp, 0, ksp->rnorm, &ksp->reason, ksp->cnvP);
 44:   if (ksp->reason) return 0;

 46:   /* Make the initial Rp == R */
 47:   VecCopy(R, RP);

 49:   /* Set the initial conditions */
 50:   etaold = 0.0;
 51:   psiold = 0.0;
 52:   tau    = dp;
 53:   dpold  = dp;

 55:   VecDot(R, RP, &rhoold); /* rhoold = (r,rp)     */
 56:   VecCopy(R, U);
 57:   VecCopy(R, P);
 58:   KSP_PCApplyBAorAB(ksp, P, V, T);
 59:   VecSet(D, 0.0);

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

 90:       dpest = PetscSqrtReal(2 * i + m + 2.0) * tau;
 91:       PetscObjectSAWsTakeAccess((PetscObject)ksp);
 92:       if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dpest;
 93:       else ksp->rnorm = 0.0;
 94:       PetscObjectSAWsGrantAccess((PetscObject)ksp);
 95:       KSPLogResidualHistory(ksp, ksp->rnorm);
 96:       KSPMonitor(ksp, i + 1, ksp->rnorm);
 97:       (*ksp->converged)(ksp, i + 1, ksp->rnorm, &ksp->reason, ksp->cnvP);
 98:       if (ksp->reason) break;

100:       etaold = eta;
101:       psiold = psi;
102:     }
103:     if (ksp->reason) break;

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

112:     rhoold = rho;
113:     dpold  = dp;

115:     i++;
116:   } while (i < ksp->max_it);
117:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

119:   KSPUnwindPreconditioner(ksp, X, T);
120:   return 0;
121: }

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

126:    Level: beginner

128:    Notes:
129:    Supports left and right preconditioning, but not symmetric

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

135:    References:
136: .  * - Freund, 1993

138: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPTCQMR`
139: M*/
140: PETSC_EXTERN PetscErrorCode KSPCreate_TFQMR(KSP ksp)
141: {
142:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
143:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
144:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);
145:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 1);

147:   ksp->data                = (void *)0;
148:   ksp->ops->setup          = KSPSetUp_TFQMR;
149:   ksp->ops->solve          = KSPSolve_TFQMR;
150:   ksp->ops->destroy        = KSPDestroyDefault;
151:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
152:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
153:   ksp->ops->setfromoptions = NULL;
154:   ksp->ops->view           = NULL;
155:   return 0;
156: }