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: if (ksp->pc_side == PC_RIGHT) {
15: SETERRQ(2,"No right preconditioning for KSPLCD");
16: } else if (ksp->pc_side == PC_SYMMETRIC) {
17: SETERRQ(2,"No symmetric preconditioning for KSPLCD");
18: }
20: /* get work vectors needed by LCD */
21: KSPDefaultGetWork(ksp,2);
22:
23: VecDuplicateVecs(ksp->work[0],restart+1,&lcd->P);
24: VecDuplicateVecs(ksp->work[0], restart + 1, &lcd->Q);
25: PetscLogObjectMemory(ksp,2*(restart+2)*sizeof(Vec));
26: return(0);
27: }
29: /* KSPSolve_LCD - This routine actually applies the left conjugate
30: direction method
32: Input Parameter:
33: . ksp - the Krylov space object that was set to use LCD, by, for
34: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPLCD);
36: Output Parameter:
37: . its - number of iterations used
39: */
42: PetscErrorCode KSPSolve_LCD(KSP ksp)
43: {
45: PetscInt it,j,max_k;
46: PetscScalar alfa, beta, num, den, mone, pone;
47: PetscReal rnorm;
48: Vec X,B,R,Z;
49: KSP_LCD *lcd;
50: Mat Amat,Pmat;
51: MatStructure pflag;
52: PetscTruth diagonalscale;
55:
56: PCDiagonalScale(ksp->pc,&diagonalscale);
57: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
59: lcd = (KSP_LCD*)ksp->data;
60: X = ksp->vec_sol;
61: B = ksp->vec_rhs;
62: R = ksp->work[0];
63: Z = ksp->work[1];
64: max_k = lcd->restart;
65: mone = -1;
66: pone = 1;
68: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
70: ksp->its = 0;
71: if (!ksp->guess_zero) {
72: KSP_MatMult(ksp,Amat,X,Z); /* z <- b - Ax */
73: VecAYPX(Z,mone,B);
74: } else {
75: VecCopy(B,Z); /* z <- b (x is 0) */
76: }
77:
78: KSP_PCApply(ksp,Z,R); /* r <- M^-1z */
79: VecNorm(R,NORM_2,&rnorm);
80: KSPLogResidualHistory(ksp,rnorm);
81: KSPMonitor(ksp,0,rnorm); /* call any registered monitor routines */
82: ksp->rnorm = rnorm;
83:
84: /* test for convergence */
85: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
86: if (ksp->reason) return(0);
88: it = 0;
89: VecCopy(R,lcd->P[0]);
90:
91: while (!ksp->reason && ksp->its < ksp->max_it) {
92: it = 0;
93: KSP_MatMult(ksp,Amat,lcd->P[it],Z);
94: KSP_PCApply(ksp,Z,lcd->Q[it]);
95:
96: while(!ksp->reason && it < max_k && ksp->its < ksp->max_it) {
97: ksp->its++;
98: VecDot(lcd->P[it],R,&num);
99: VecDot(lcd->P[it],lcd->Q[it], &den);
100: alfa = num/den;
101: VecAXPY(X,alfa,lcd->P[it]);
102: VecAXPY(R,-alfa,lcd->Q[it]);
103: VecNorm(R,NORM_2,&rnorm);
105: ksp->rnorm = rnorm;
106: KSPLogResidualHistory(ksp,rnorm);
107: KSPMonitor(ksp,ksp->its,rnorm);
108: (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
109:
110: if (ksp->reason) break;
111:
112: VecCopy(R,lcd->P[it+1]);
113: KSP_MatMult(ksp,Amat,lcd->P[it+1],Z);
114: KSP_PCApply(ksp,Z,lcd->Q[it+1]);
115:
116: for( j = 0; j <= it; j++) {
117: VecDot(lcd->P[j],lcd->Q[it+1],&num);
118: VecDot(lcd->P[j],lcd->Q[j],&den);
119: beta = - num/den;
120: VecAXPY(lcd->P[it+1],beta,lcd->P[j]);
121: VecAXPY(lcd->Q[it+1],beta,lcd->Q[j]);
122: }
123: it++;
124: }
125: VecCopy(lcd->P[it],lcd->P[0]);
126: }
127: if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
128: VecCopy(X,ksp->vec_sol);
129:
130: return(0);
131: }
132: /*
133: KSPDestroy_LCD - Frees all memory space used by the Krylov method
135: */
138: PetscErrorCode KSPDestroy_LCD(KSP ksp)
139: {
140: KSP_LCD *lcd = (KSP_LCD*)ksp->data;
144: KSPDefaultFreeWork(ksp);
145: if (lcd->P) { VecDestroyVecs(lcd->P, lcd->restart+1); }
146: if (lcd->Q) { VecDestroyVecs(lcd->Q, lcd->restart+1); }
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: PetscTruth iascii;
169: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&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: } else {
174: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP LCD",((PetscObject)viewer)->type_name);
175: }
176: return(0);
177: }
179: /*
180: KSPSetFromOptions_LCD - Checks the options database for options related to the
181: LCD method.
182: */
185: PetscErrorCode KSPSetFromOptions_LCD(KSP ksp)
186: {
188: PetscTruth flg;
189: KSP_LCD *lcd = (KSP_LCD *)ksp->data;
190:
192: PetscOptionsHead("KSP LCD options");
193: PetscOptionsInt("-ksp_lcd_restart","Number of vectors conjugate","KSPLCDSetRestart",lcd->restart,&lcd->restart,&flg);
194: if(flg && lcd->restart < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
195: PetscOptionsReal("-ksp_lcd_haptol","Tolerance for exact convergence (happy ending)","KSPLCDSetHapTol",lcd->haptol,&lcd->haptol,&flg);
196: if (flg && lcd->haptol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
197: return(0);
198: }
200: /*MC
201: KSPLCD - Implements the LCD (left conjugate direction) method in PETSc.
203: Options Database Keys:
204: + -ksp_lcd_restart - number of vectors conjudate
205: - -ksp_lcd_haptol - tolerance for exact convergence (happing ending)
207: Level: beginner
210: References:
211: - J.Y. Yuan, G.H.Golub, R.J. Plemmons, and W.A.G. Cecilio. Semiconjugate
212: direction methods for real positive definite system. BIT Numerical
213: Mathematics, 44(1):189-207,2004.
214: - Y. Dai and J.Y. Yuan. Study on semi-conjugate direction methods for
215: non-symmetric systems. International Journal for Numerical Methods in
216: Engineering, 60:1383-1399,2004.
217: - L. Catabriga, A.L.G.A. Coutinho, and L.P.Franca. Evaluating the LCD
218: algorithm for solving linear systems of equations arising from implicit
219: SUPG formulation of compressible flows. International Journal for
220: Numerical Methods in Engineering, 60:1513-1534,2004
221: - L. Catabriga, A. M. P. Valli, B. Z. Melotti, L. M. Pessoa,
222: A. L. G. A. Coutinho, Performance of LCD iterative method in the finite
223: element and finite difference solution of convection-diffusion
224: equations, Communications in Numerical Methods in Engineering, (Early
225: View).
227: Contributed by: Lucia Catabriga <luciac@ices.utexas.edu>
230: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
231: KSPCGSetType(), KSPLCDSetRestart(), KSPLCDSetHapTol()
233: M*/
238: PetscErrorCode KSPCreate_LCD(KSP ksp)
239: {
241: KSP_LCD *lcd;
244: PetscNewLog(ksp,KSP_LCD,&lcd);
245: ksp->data = (void*)lcd;
246: ksp->pc_side = PC_LEFT;
247: lcd->restart = 30;
248: lcd->haptol = 1.0e-30;
250: /*
251: Sets the functions that are associated with this data structure
252: (in C++ this is the same as defining virtual functions)
253: */
254: ksp->ops->setup = KSPSetUp_LCD;
255: ksp->ops->solve = KSPSolve_LCD;
256: ksp->ops->destroy = KSPDestroy_LCD;
257: ksp->ops->view = KSPView_LCD;
258: ksp->ops->setfromoptions = KSPSetFromOptions_LCD;
259: ksp->ops->buildsolution = KSPDefaultBuildSolution;
260: ksp->ops->buildresidual = KSPDefaultBuildResidual;
262: return(0);
263: }