Actual source code: lcd.c

  1: #define PETSCKSP_DLL

 3:  #include src/ksp/ksp/impls/lcd/lcdctx.h

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

 14:   /* 
 15:        This implementation of LCD only handles left preconditioning
 16:      so generate an error otherwise.
 17:   */
 18:   if (ksp->pc_side == PC_RIGHT) {
 19:     SETERRQ(2,"No right preconditioning for KSPLCD");
 20:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 21:     SETERRQ(2,"No symmetric preconditioning for KSPLCD");
 22:   }

 24:   /* get work vectors needed by LCD */
 25:   KSPDefaultGetWork(ksp,2);
 26: 
 27:   VecDuplicateVecs(ksp->work[0],restart+1,&lcd->P);
 28:   VecDuplicateVecs(ksp->work[0], restart + 1, &lcd->Q);
 29:   PetscLogObjectMemory(ksp,2*(restart+2)*sizeof(Vec));
 30:   return(0);
 31: }

 33: /*     KSPSolve_LCD - This routine actually applies the left conjugate
 34:     direction method

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

 40:    Output Parameter:
 41: .     its - number of iterations used

 43: */
 46: PetscErrorCode  KSPSolve_LCD(KSP ksp)
 47: {
 49:   PetscInt       it,j,max_k;
 50:   PetscScalar    alfa, beta, num, den, mone, pone;
 51:   PetscReal      rnorm;
 52:   Vec            X,B,R,Z;
 53:   KSP_LCD        *lcd;
 54:   Mat            Amat,Pmat;
 55:   MatStructure   pflag;
 56:   PetscTruth     diagonalscale;

 59: 
 60:   PCDiagonalScale(ksp->pc,&diagonalscale);
 61:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 63:   lcd            = (KSP_LCD*)ksp->data;
 64:   X              = ksp->vec_sol;
 65:   B              = ksp->vec_rhs;
 66:   R              = ksp->work[0];
 67:   Z              = ksp->work[1];
 68:   max_k          = lcd->restart;
 69:   mone = -1;
 70:   pone = 1;

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

 74:   ksp->its = 0;
 75:   if (!ksp->guess_zero) {
 76:     KSP_MatMult(ksp,Amat,X,Z);             /*   z <- b - Ax       */
 77:     VecAYPX(Z,mone,B);
 78:   } else {
 79:     VecCopy(B,Z);                         /*     z <- b (x is 0) */
 80:   }
 81: 
 82:   KSP_PCApply(ksp,Z,R);                   /*     r <- M^-1z         */
 83:   VecNorm(R,NORM_2,&rnorm);
 84:   KSPLogResidualHistory(ksp,rnorm);
 85:   KSPMonitor(ksp,0,rnorm);                    /* call any registered monitor routines */
 86:   ksp->rnorm = rnorm;
 87: 
 88:    /* test for convergence */
 89:   (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
 90:   if (ksp->reason) return(0);

 92:   it = 0;
 93:   VecCopy(R,lcd->P[0]);
 94: 
 95:   while (!ksp->reason && ksp->its < ksp->max_it) {
 96:     it = 0;
 97:     KSP_MatMult(ksp,Amat,lcd->P[it],Z);
 98:     KSP_PCApply(ksp,Z,lcd->Q[it]);
 99: 
100:     while(!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
101:       ksp->its++;
102:       VecDot(lcd->P[it],R,&num);
103:       VecDot(lcd->P[it],lcd->Q[it], &den);
104:       alfa = num/den;
105:       VecAXPY(X,alfa,lcd->P[it]);
106:       VecAXPY(R,-alfa,lcd->Q[it]);
107:       VecNorm(R,NORM_2,&rnorm);

109:       ksp->rnorm = rnorm;
110:       KSPLogResidualHistory(ksp,rnorm);
111:       KSPMonitor(ksp,ksp->its,rnorm);
112:       (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
113: 
114:       if (ksp->reason) break;
115: 
116:       VecCopy(R,lcd->P[it+1]);
117:       KSP_MatMult(ksp,Amat,lcd->P[it+1],Z);
118:       KSP_PCApply(ksp,Z,lcd->Q[it+1]);
119: 
120:       for( j = 0; j <= it; j++)        {
121:         VecDot(lcd->P[j],lcd->Q[it+1],&num);
122:         VecDot(lcd->P[j],lcd->Q[j],&den);
123:         beta = - num/den;
124:         VecAXPY(lcd->P[it+1],beta,lcd->P[j]);
125:         VecAXPY(lcd->Q[it+1],beta,lcd->Q[j]);
126:       }
127:       it++;
128:     }
129:     VecCopy(lcd->P[it],lcd->P[0]);
130:   }
131:   if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
132:   VecCopy(X,ksp->vec_sol);
133: 
134:   return(0);
135: }
136: /*
137:        KSPDestroy_LCD - Frees all memory space used by the Krylov method

139: */
142: PetscErrorCode KSPDestroy_LCD(KSP ksp)
143: {
144:   KSP_LCD         *lcd = (KSP_LCD*)ksp->data;

148:   KSPDefaultFreeWork(ksp);
149:   if (lcd->P) { VecDestroyVecs(lcd->P, lcd->restart+1); }
150:   if (lcd->Q) { VecDestroyVecs(lcd->Q, lcd->restart+1); }
151:   PetscFree(ksp->data);
152:   return(0);
153: }

155: /*
156:      KSPView_LCD - Prints information about the current Krylov method being used

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

162: */
165: PetscErrorCode KSPView_LCD(KSP ksp,PetscViewer viewer)
166: {

168:   KSP_LCD         *lcd = (KSP_LCD *)ksp->data;
170:   PetscTruth     iascii;

173:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
174:   if (iascii) {
175:       PetscViewerASCIIPrintf(viewer,"  LCD: restart=%d\n",lcd->restart);
176:       PetscViewerASCIIPrintf(viewer,"  LCD: happy breakdown tolerance %g\n",lcd->haptol);
177:   } else {
178:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP LCD",((PetscObject)viewer)->type_name);
179:   }
180:   return(0);
181: }

183: /*
184:     KSPSetFromOptions_LCD - Checks the options database for options related to the 
185:                             LCD method.
186: */
189: PetscErrorCode KSPSetFromOptions_LCD(KSP ksp)
190: {
192:   PetscTruth     flg;
193:   KSP_LCD        *lcd = (KSP_LCD *)ksp->data;
194: 
196:   PetscOptionsHead("KSP LCD options");
197:   PetscOptionsInt("-ksp_lcd_restart","Number of vectors conjugate","KSPLCDSetRestart",lcd->restart,&lcd->restart,&flg);
198:   if(flg && lcd->restart < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
199:   PetscOptionsReal("-ksp_lcd_haptol","Tolerance for exact convergence (happy ending)","KSPLCDSetHapTol",lcd->haptol,&lcd->haptol,&flg);
200:   if (flg && lcd->haptol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
201:   return(0);
202: }

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

207:    Options Database Keys:
208: +   -ksp_lcd_restart - number of vectors conjudate
209: -   -ksp_lcd_haptol - tolerance for exact convergence (happing ending)

211:    Level: beginner


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

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


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

237: M*/

242: PetscErrorCode KSPCreate_LCD(KSP ksp)
243: {
245:   KSP_LCD         *lcd;

248:   PetscNewLog(ksp,KSP_LCD,&lcd);
249:   ksp->data                      = (void*)lcd;
250:   ksp->pc_side                   = PC_LEFT;
251:   lcd->restart                   = 30;
252:   lcd->haptol                    = 1.0e-30;

254:   /*
255:        Sets the functions that are associated with this data structure 
256:        (in C++ this is the same as defining virtual functions)
257:   */
258:   ksp->ops->setup                = KSPSetUp_LCD;
259:   ksp->ops->solve                = KSPSolve_LCD;
260:   ksp->ops->destroy              = KSPDestroy_LCD;
261:   ksp->ops->view                 = KSPView_LCD;
262:   ksp->ops->setfromoptions       = KSPSetFromOptions_LCD;
263:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
264:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

266:   return(0);
267: }