Actual source code: lcd.c

petsc-3.5.1 2014-08-06
Report Typos and Errors
  2: #include <../src/ksp/ksp/impls/lcd/lcdimpl.h>

  6: PetscErrorCode KSPSetUp_LCD(KSP ksp)
  7: {
  8:   KSP_LCD        *lcd = (KSP_LCD*)ksp->data;
 10:   PetscInt       restart = lcd->restart;

 13:   /* get work vectors needed by LCD */
 14:   KSPSetWorkVecs(ksp,2);

 16:   VecDuplicateVecs(ksp->work[0],restart+1,&lcd->P);
 17:   VecDuplicateVecs(ksp->work[0], restart + 1, &lcd->Q);
 18:   PetscLogObjectMemory((PetscObject)ksp,2*(restart+2)*sizeof(Vec));
 19:   return(0);
 20: }

 22: /*     KSPSolve_LCD - This routine actually applies the left conjugate
 23:     direction method

 25:    Input Parameter:
 26: .     ksp - the Krylov space object that was set to use LCD, by, for
 27:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPLCD);

 29:    Output Parameter:
 30: .     its - number of iterations used

 32: */
 35: PetscErrorCode  KSPSolve_LCD(KSP ksp)
 36: {
 38:   PetscInt       it,j,max_k;
 39:   PetscScalar    alfa, beta, num, den, mone;
 40:   PetscReal      rnorm;
 41:   Vec            X,B,R,Z;
 42:   KSP_LCD        *lcd;
 43:   Mat            Amat,Pmat;
 44:   PetscBool      diagonalscale;

 47:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 48:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 50:   lcd   = (KSP_LCD*)ksp->data;
 51:   X     = ksp->vec_sol;
 52:   B     = ksp->vec_rhs;
 53:   R     = ksp->work[0];
 54:   Z     = ksp->work[1];
 55:   max_k = lcd->restart;
 56:   mone  = -1;

 58:   PCGetOperators(ksp->pc,&Amat,&Pmat);

 60:   ksp->its = 0;
 61:   if (!ksp->guess_zero) {
 62:     KSP_MatMult(ksp,Amat,X,Z);             /*   z <- b - Ax       */
 63:     VecAYPX(Z,mone,B);
 64:   } else {
 65:     VecCopy(B,Z);                         /*     z <- b (x is 0) */
 66:   }

 68:   KSP_PCApply(ksp,Z,R);                   /*     r <- M^-1z         */
 69:   VecNorm(R,NORM_2,&rnorm);
 70:   KSPLogResidualHistory(ksp,rnorm);
 71:   KSPMonitor(ksp,0,rnorm);
 72:   ksp->rnorm = rnorm;

 74:   /* test for convergence */
 75:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
 76:   if (ksp->reason) return(0);

 78:   it = 0;
 79:   VecCopy(R,lcd->P[0]);

 81:   while (!ksp->reason && ksp->its < ksp->max_it) {
 82:     it   = 0;
 83:     KSP_MatMult(ksp,Amat,lcd->P[it],Z);
 84:     KSP_PCApply(ksp,Z,lcd->Q[it]);

 86:     while (!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
 87:       ksp->its++;
 88:       VecDot(lcd->P[it],R,&num);
 89:       VecDot(lcd->P[it],lcd->Q[it], &den);
 90:       alfa = num/den;
 91:       VecAXPY(X,alfa,lcd->P[it]);
 92:       VecAXPY(R,-alfa,lcd->Q[it]);
 93:       VecNorm(R,NORM_2,&rnorm);

 95:       ksp->rnorm = rnorm;
 96:       KSPLogResidualHistory(ksp,rnorm);
 97:       KSPMonitor(ksp,ksp->its,rnorm);
 98:       (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);

100:       if (ksp->reason) break;

102:       VecCopy(R,lcd->P[it+1]);
103:       KSP_MatMult(ksp,Amat,lcd->P[it+1],Z);
104:       KSP_PCApply(ksp,Z,lcd->Q[it+1]);

106:       for (j = 0; j <= it; j++) {
107:         VecDot(lcd->P[j],lcd->Q[it+1],&num);
108:         VecDot(lcd->P[j],lcd->Q[j],&den);
109:         beta = -num/den;
110:         VecAXPY(lcd->P[it+1],beta,lcd->P[j]);
111:         VecAXPY(lcd->Q[it+1],beta,lcd->Q[j]);
112:       }
113:       it++;
114:     }
115:     VecCopy(lcd->P[it],lcd->P[0]);
116:   }
117:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
118:   VecCopy(X,ksp->vec_sol);
119:   return(0);
120: }
121: /*
122:        KSPDestroy_LCD - Frees all memory space used by the Krylov method

124: */
127: PetscErrorCode KSPReset_LCD(KSP ksp)
128: {
129:   KSP_LCD        *lcd = (KSP_LCD*)ksp->data;

133:   if (lcd->P) { VecDestroyVecs(lcd->restart+1,&lcd->P);}
134:   if (lcd->Q) { VecDestroyVecs(lcd->restart+1,&lcd->Q);}
135:   return(0);
136: }


141: PetscErrorCode KSPDestroy_LCD(KSP ksp)
142: {

146:   KSPReset_LCD(ksp);
147:   PetscFree(ksp->data);
148:   return(0);
149: }

151: /*
152:      KSPView_LCD - Prints information about the current Krylov method being used

154:       Currently this only prints information to a file (or stdout) about the
155:       symmetry of the problem. If your Krylov method has special options or
156:       flags that information should be printed here.

158: */
161: PetscErrorCode KSPView_LCD(KSP ksp,PetscViewer viewer)
162: {

164:   KSP_LCD        *lcd = (KSP_LCD*)ksp->data;
166:   PetscBool      iascii;

169:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
170:   if (iascii) {
171:     PetscViewerASCIIPrintf(viewer,"  LCD: restart=%d\n",lcd->restart);
172:     PetscViewerASCIIPrintf(viewer,"  LCD: happy breakdown tolerance %g\n",lcd->haptol);
173:   }
174:   return(0);
175: }

177: /*
178:     KSPSetFromOptions_LCD - Checks the options database for options related to the
179:                             LCD method.
180: */
183: PetscErrorCode KSPSetFromOptions_LCD(KSP ksp)
184: {
186:   PetscBool      flg;
187:   KSP_LCD        *lcd = (KSP_LCD*)ksp->data;

190:   PetscOptionsHead("KSP LCD options");
191:   PetscOptionsInt("-ksp_lcd_restart","Number of vectors conjugate","KSPLCDSetRestart",lcd->restart,&lcd->restart,&flg);
192:   if (flg && lcd->restart < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
193:   PetscOptionsReal("-ksp_lcd_haptol","Tolerance for exact convergence (happy ending)","KSPLCDSetHapTol",lcd->haptol,&lcd->haptol,&flg);
194:   if (flg && lcd->haptol < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
195:   return(0);
196: }

198: /*MC
199:      KSPLCD -  Implements the LCD (left conjugate direction) method in PETSc.

201:    Options Database Keys:
202: +   -ksp_lcd_restart - number of vectors conjudate
203: -   -ksp_lcd_haptol - tolerance for exact convergence (happing ending)

205:    Level: beginner

207:     Notes: Support only for left preconditioning

209:     References:
210:    - J.Y. Yuan, G.H.Golub, R.J. Plemmons, and W.A.G. Cecilio. Semiconjugate
211:      direction methods for real positive definite system. BIT Numerical
212:      Mathematics, 44(1):189-207,2004.
213:    - Y. Dai and J.Y. Yuan. Study on semi-conjugate direction methods for
214:      non-symmetric systems. International Journal for Numerical Methods in
215:      Engineering, 60:1383-1399,2004.
216:    - L. Catabriga, A.L.G.A. Coutinho, and L.P.Franca. Evaluating the LCD
217:      algorithm for solving linear systems of equations arising from implicit
218:      SUPG formulation of compressible flows. International Journal for
219:      Numerical Methods in Engineering, 60:1513-1534,2004
220:    - L. Catabriga, A. M. P. Valli, B. Z. Melotti, L. M. Pessoa,
221:      A. L. G. A. Coutinho, Performance of LCD iterative method in the finite
222:      element and finite difference solution of convection-diffusion
223:      equations,  Communications in Numerical Methods in Engineering, (Early
224:      View).

226:   Contributed by: Lucia Catabriga <luciac@ices.utexas.edu>


229: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
230:            KSPCGSetType(), KSPLCDSetRestart(), KSPLCDSetHapTol()

232: M*/

236: PETSC_EXTERN PetscErrorCode KSPCreate_LCD(KSP ksp)
237: {
239:   KSP_LCD        *lcd;

242:   PetscNewLog(ksp,&lcd);
243:   ksp->data    = (void*)lcd;
244:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
245:   lcd->restart = 30;
246:   lcd->haptol  = 1.0e-30;

248:   /*
249:        Sets the functions that are associated with this data structure
250:        (in C++ this is the same as defining virtual functions)
251:   */
252:   ksp->ops->setup          = KSPSetUp_LCD;
253:   ksp->ops->solve          = KSPSolve_LCD;
254:   ksp->ops->reset          = KSPReset_LCD;
255:   ksp->ops->destroy        = KSPDestroy_LCD;
256:   ksp->ops->view           = KSPView_LCD;
257:   ksp->ops->setfromoptions = KSPSetFromOptions_LCD;
258:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
259:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
260:   return(0);
261: }