Actual source code: stcg.c
1: #define PETSCKSP_DLL
3: #include private/kspimpl.h
4: #include ../src/ksp/ksp/impls/cg/stcg/stcgimpl.h
6: #define STCG_PRECONDITIONED_DIRECTION 0
7: #define STCG_UNPRECONDITIONED_DIRECTION 1
8: #define STCG_DIRECTION_TYPES 2
10: static const char *DType_Table[64] = {"preconditioned", "unpreconditioned"};
14: /*@
15: KSPSTCGSetRadius - Sets the radius of the trust region.
17: Collective on KSP
19: Input Parameters:
20: + ksp - the iterative context
21: - radius - the trust region radius (Infinity is the default)
23: Options Database Key:
24: . -ksp_stcg_radius <r>
26: Level: advanced
28: .keywords: KSP, STCG, set, trust region radius
29: @*/
30: PetscErrorCode KSPSTCGSetRadius(KSP ksp, PetscReal radius)
31: {
32: PetscErrorCode ierr, (*f)(KSP, PetscReal);
36: if (radius < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Radius negative");
37: PetscObjectQueryFunction((PetscObject)ksp, "KSPSTCGSetRadius_C", (void (**)(void))&f);
38: if (f) {
39: (*f)(ksp, radius);
40: }
41: return(0);
42: }
46: /*@
47: KSPSTCGGetNormD - Got norm of the direction.
49: Collective on KSP
51: Input Parameters:
52: + ksp - the iterative context
53: - norm_d - the norm of the direction
55: Level: advanced
57: .keywords: KSP, STCG, get, norm direction
58: @*/
59: PetscErrorCode KSPSTCGGetNormD(KSP ksp, PetscReal *norm_d)
60: {
61: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
65: PetscObjectQueryFunction((PetscObject)ksp, "KSPSTCGGetNormD_C", (void (**)(void))&f);
66: if (f) {
67: (*f)(ksp, norm_d);
68: }
69: return(0);
70: }
74: /*@
75: KSPSTCGGetObjFcn - Get objective function value.
77: Collective on KSP
79: Input Parameters:
80: + ksp - the iterative context
81: - o_fcn - the objective function value
83: Level: advanced
85: .keywords: KSP, STCG, get, objective function
86: @*/
87: PetscErrorCode KSPSTCGGetObjFcn(KSP ksp, PetscReal *o_fcn)
88: {
89: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
93: PetscObjectQueryFunction((PetscObject)ksp, "KSPSTCGGetObjFcn_C", (void (**)(void))&f);
94: if (f) {
95: (*f)(ksp, o_fcn);
96: }
97: return(0);
98: }
102: PetscErrorCode KSPSolve_STCG(KSP ksp)
103: {
104: #ifdef PETSC_USE_COMPLEX
105: SETERRQ(PETSC_ERR_SUP, "STCG is not available for complex systems");
106: #else
107: KSP_STCG *cg = (KSP_STCG *)ksp->data;
110: MatStructure pflag;
111: Mat Qmat, Mmat;
112: Vec r, z, p, d;
113: PC pc;
115: PetscReal norm_r, norm_d, norm_dp1, norm_p, dMp;
116: PetscReal alpha, beta, kappa, rz, rzm1;
117: PetscReal rr, r2, step;
119: PetscInt max_cg_its;
121: PetscTruth diagonalscale;
124: /***************************************************************************/
125: /* Check the arguments and parameters. */
126: /***************************************************************************/
128: PCDiagonalScale(ksp->pc, &diagonalscale);
129: if (diagonalscale) {
130: SETERRQ1(PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
131: }
133: if (cg->radius < 0.0) {
134: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Input error: radius < 0");
135: }
137: /***************************************************************************/
138: /* Get the workspace vectors and initialize variables */
139: /***************************************************************************/
141: r2 = cg->radius * cg->radius;
142: r = ksp->work[0];
143: z = ksp->work[1];
144: p = ksp->work[2];
145: d = ksp->vec_sol;
146: pc = ksp->pc;
148: PCGetOperators(pc, &Qmat, &Mmat, &pflag);
150: VecGetSize(d, &max_cg_its);
151: max_cg_its = PetscMin(max_cg_its, ksp->max_it);
152: ksp->its = 0;
154: /***************************************************************************/
155: /* Initialize objective function and direction. */
156: /***************************************************************************/
158: cg->o_fcn = 0.0;
160: VecSet(d, 0.0); /* d = 0 */
161: cg->norm_d = 0.0;
163: /***************************************************************************/
164: /* Begin the conjugate gradient method. Check the right-hand side for */
165: /* numerical problems. The check for not-a-number and infinite values */
166: /* need be performed only once. */
167: /***************************************************************************/
169: VecCopy(ksp->vec_rhs, r); /* r = -grad */
170: VecDot(r, r, &rr); /* rr = r^T r */
171: if (PetscIsInfOrNanScalar(rr)) {
172: /*************************************************************************/
173: /* The right-hand side contains not-a-number or an infinite value. */
174: /* The gradient step does not work; return a zero value for the step. */
175: /*************************************************************************/
177: ksp->reason = KSP_DIVERGED_NAN;
178: PetscInfo1(ksp, "KSPSolve_STCG: bad right-hand side: rr=%g\n", rr);
179: return(0);
180: }
182: /***************************************************************************/
183: /* Check the preconditioner for numerical problems and for positive */
184: /* definiteness. The check for not-a-number and infinite values need be */
185: /* performed only once. */
186: /***************************************************************************/
188: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
189: VecDot(r, z, &rz); /* rz = r^T inv(M) r */
190: if (PetscIsInfOrNanScalar(rz)) {
191: /*************************************************************************/
192: /* The preconditioner contains not-a-number or an infinite value. */
193: /* Return the gradient direction intersected with the trust region. */
194: /*************************************************************************/
196: ksp->reason = KSP_DIVERGED_NAN;
197: PetscInfo1(ksp, "KSPSolve_STCG: bad preconditioner: rz=%g\n", rz);
199: if (cg->radius != 0) {
200: if (r2 >= rr) {
201: alpha = 1.0;
202: cg->norm_d = sqrt(rr);
203: }
204: else {
205: alpha = sqrt(r2 / rr);
206: cg->norm_d = cg->radius;
207: }
209: VecAXPY(d, alpha, r); /* d = d + alpha r */
211: /***********************************************************************/
212: /* Compute objective function. */
213: /***********************************************************************/
215: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
216: VecAYPX(z, -0.5, ksp->vec_rhs);
217: VecDot(d, z, &cg->o_fcn);
218: cg->o_fcn = -cg->o_fcn;
219: ++ksp->its;
220: }
221: return(0);
222: }
224: if (rz < 0.0) {
225: /*************************************************************************/
226: /* The preconditioner is indefinite. Because this is the first */
227: /* and we do not have a direction yet, we use the gradient step. Note */
228: /* that we cannot use the preconditioned norm when computing the step */
229: /* because the matrix is indefinite. */
230: /*************************************************************************/
232: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
233: PetscInfo1(ksp, "KSPSolve_STCG: indefinite preconditioner: rz=%g\n", rz);
235: if (cg->radius != 0.0) {
236: if (r2 >= rr) {
237: alpha = 1.0;
238: cg->norm_d = sqrt(rr);
239: }
240: else {
241: alpha = sqrt(r2 / rr);
242: cg->norm_d = cg->radius;
243: }
245: VecAXPY(d, alpha, r); /* d = d + alpha r */
247: /***********************************************************************/
248: /* Compute objective function. */
249: /***********************************************************************/
251: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
252: VecAYPX(z, -0.5, ksp->vec_rhs);
253: VecDot(d, z, &cg->o_fcn);
254: cg->o_fcn = -cg->o_fcn;
255: ++ksp->its;
256: }
257: return(0);
258: }
260: /***************************************************************************/
261: /* As far as we know, the preconditioner is positive semidefinite. */
262: /* Compute and log the residual. Check convergence because this */
263: /* initializes things, but do not terminate until at least one conjugate */
264: /* gradient iteration has been performed. */
265: /***************************************************************************/
267: switch(ksp->normtype) {
268: case KSP_NORM_PRECONDITIONED:
269: VecNorm(z, NORM_2, &norm_r); /* norm_r = |z| */
270: break;
272: case KSP_NORM_UNPRECONDITIONED:
273: norm_r = sqrt(rr); /* norm_r = |r| */
274: break;
276: case KSP_NORM_NATURAL:
277: norm_r = sqrt(rz); /* norm_r = |r|_M */
278: break;
280: default:
281: norm_r = 0.0;
282: break;
283: }
285: KSPLogResidualHistory(ksp, norm_r);
286: KSPMonitor(ksp, ksp->its, norm_r);
287: ksp->rnorm = norm_r;
289: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
291: /***************************************************************************/
292: /* Compute the first direction and update the iteration. */
293: /***************************************************************************/
295: VecCopy(z, p); /* p = z */
296: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
297: ++ksp->its;
299: /***************************************************************************/
300: /* Check the matrix for numerical problems. */
301: /***************************************************************************/
303: VecDot(p, z, &kappa); /* kappa = p^T Q p */
304: if (PetscIsInfOrNanScalar(kappa)) {
305: /*************************************************************************/
306: /* The matrix produced not-a-number or an infinite value. In this case, */
307: /* we must stop and use the gradient direction. This condition need */
308: /* only be checked once. */
309: /*************************************************************************/
311: ksp->reason = KSP_DIVERGED_NAN;
312: PetscInfo1(ksp, "KSPSolve_STCG: bad matrix: kappa=%g\n", kappa);
314: if (cg->radius) {
315: if (r2 >= rr) {
316: alpha = 1.0;
317: cg->norm_d = sqrt(rr);
318: }
319: else {
320: alpha = sqrt(r2 / rr);
321: cg->norm_d = cg->radius;
322: }
324: VecAXPY(d, alpha, r); /* d = d + alpha r */
326: /***********************************************************************/
327: /* Compute objective function. */
328: /***********************************************************************/
330: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
331: VecAYPX(z, -0.5, ksp->vec_rhs);
332: VecDot(d, z, &cg->o_fcn);
333: cg->o_fcn = -cg->o_fcn;
334: ++ksp->its;
335: }
336: return(0);
337: }
339: /***************************************************************************/
340: /* Initialize variables for calculating the norm of the direction. */
341: /***************************************************************************/
343: dMp = 0.0;
344: norm_d = 0.0;
345: switch(cg->dtype) {
346: case STCG_PRECONDITIONED_DIRECTION:
347: norm_p = rz;
348: break;
350: default:
351: VecDot(p, p, &norm_p);
352: break;
353: }
355: /***************************************************************************/
356: /* Check for negative curvature. */
357: /***************************************************************************/
359: if (kappa <= 0.0) {
360: /*************************************************************************/
361: /* In this case, the matrix is indefinite and we have encountered a */
362: /* direction of negative curvature. Because negative curvature occurs */
363: /* during the first step, we must follow a direction. */
364: /*************************************************************************/
366: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
367: PetscInfo1(ksp, "KSPSolve_STCG: negative curvature: kappa=%g\n", kappa);
369: if (cg->radius != 0.0 && norm_p > 0.0) {
370: /***********************************************************************/
371: /* Follow direction of negative curvature to the boundary of the */
372: /* trust region. */
373: /***********************************************************************/
375: step = sqrt(r2 / norm_p);
376: cg->norm_d = cg->radius;
378: VecAXPY(d, step, p); /* d = d + step p */
380: /***********************************************************************/
381: /* Update objective function. */
382: /***********************************************************************/
384: cg->o_fcn += step * (0.5 * step * kappa - rz);
385: }
386: else if (cg->radius != 0.0) {
387: /***********************************************************************/
388: /* The norm of the preconditioned direction is zero; use the gradient */
389: /* step. */
390: /***********************************************************************/
392: if (r2 >= rr) {
393: alpha = 1.0;
394: cg->norm_d = sqrt(rr);
395: }
396: else {
397: alpha = sqrt(r2 / rr);
398: cg->norm_d = cg->radius;
399: }
401: VecAXPY(d, alpha, r); /* d = d + alpha r */
403: /***********************************************************************/
404: /* Compute objective function. */
405: /***********************************************************************/
407: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
408: VecAYPX(z, -0.5, ksp->vec_rhs);
409: VecDot(d, z, &cg->o_fcn);
410: cg->o_fcn = -cg->o_fcn;
411: ++ksp->its;
412: }
413: return(0);
414: }
416: /***************************************************************************/
417: /* Run the conjugate gradient method until either the problem is solved, */
418: /* we encounter the boundary of the trust region, or the conjugate */
419: /* gradient method breaks down. */
420: /***************************************************************************/
422: while(1) {
423: /*************************************************************************/
424: /* Know that kappa is nonzero, because we have not broken down, so we */
425: /* can compute the steplength. */
426: /*************************************************************************/
428: alpha = rz / kappa;
430: /*************************************************************************/
431: /* Compute the steplength and check for intersection with the trust */
432: /* region. */
433: /*************************************************************************/
435: norm_dp1 = norm_d + alpha*(2.0*dMp + alpha*norm_p);
436: if (cg->radius != 0.0 && norm_dp1 >= r2) {
437: /***********************************************************************/
438: /* In this case, the matrix is positive definite as far as we know. */
439: /* However, the full step goes beyond the trust region. */
440: /***********************************************************************/
442: ksp->reason = KSP_CONVERGED_CG_CONSTRAINED;
443: PetscInfo1(ksp, "KSPSolve_STCG: constrained step: radius=%g\n", cg->radius);
445: if (norm_p > 0.0) {
446: /*********************************************************************/
447: /* Follow the direction to the boundary of the trust region. */
448: /*********************************************************************/
450: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
451: cg->norm_d = cg->radius;
453: VecAXPY(d, step, p); /* d = d + step p */
455: /*********************************************************************/
456: /* Update objective function. */
457: /*********************************************************************/
459: cg->o_fcn += step * (0.5 * step * kappa - rz);
460: }
461: else {
462: /*********************************************************************/
463: /* The norm of the direction is zero; there is nothing to follow. */
464: /*********************************************************************/
465: }
466: break;
467: }
469: /*************************************************************************/
470: /* Now we can update the direction and residual. */
471: /*************************************************************************/
473: VecAXPY(d, alpha, p); /* d = d + alpha p */
474: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
475: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
477: switch(cg->dtype) {
478: case STCG_PRECONDITIONED_DIRECTION:
479: norm_d = norm_dp1;
480: break;
482: default:
483: VecDot(d, d, &norm_d);
484: break;
485: }
486: cg->norm_d = sqrt(norm_d);
488: /*************************************************************************/
489: /* Update objective function. */
490: /*************************************************************************/
492: cg->o_fcn -= 0.5 * alpha * rz;
494: /*************************************************************************/
495: /* Check that the preconditioner appears positive semidefinite. */
496: /*************************************************************************/
498: rzm1 = rz;
499: VecDot(r, z, &rz); /* rz = r^T z */
500: if (rz < 0.0) {
501: /***********************************************************************/
502: /* The preconditioner is indefinite. */
503: /***********************************************************************/
505: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
506: PetscInfo1(ksp, "KSPSolve_STCG: cg indefinite preconditioner: rz=%g\n", rz);
507: break;
508: }
510: /*************************************************************************/
511: /* As far as we know, the preconditioner is positive semidefinite. */
512: /* Compute the residual and check for convergence. */
513: /*************************************************************************/
515: switch(ksp->normtype) {
516: case KSP_NORM_PRECONDITIONED:
517: VecNorm(z, NORM_2, &norm_r);/* norm_r = |z| */
518: break;
520: case KSP_NORM_UNPRECONDITIONED:
521: VecNorm(r, NORM_2, &norm_r);/* norm_r = |r| */
522: break;
524: case KSP_NORM_NATURAL:
525: norm_r = sqrt(rz); /* norm_r = |r|_M */
526: break;
528: default:
529: norm_r = 0.0;
530: break;
531: }
533: KSPLogResidualHistory(ksp, norm_r);
534: KSPMonitor(ksp, ksp->its, norm_r);
535: ksp->rnorm = norm_r;
536:
537: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
538: if (ksp->reason) {
539: /***********************************************************************/
540: /* The method has converged. */
541: /***********************************************************************/
543: PetscInfo2(ksp, "KSPSolve_STCG: truncated step: rnorm=%g, radius=%g\n", norm_r, cg->radius);
544: break;
545: }
547: /*************************************************************************/
548: /* We have not converged yet. Check for breakdown. */
549: /*************************************************************************/
551: beta = rz / rzm1;
552: if (fabs(beta) <= 0.0) {
553: /***********************************************************************/
554: /* Conjugate gradients has broken down. */
555: /***********************************************************************/
557: ksp->reason = KSP_DIVERGED_BREAKDOWN;
558: PetscInfo1(ksp, "KSPSolve_STCG: breakdown: beta=%g\n", beta);
559: break;
560: }
562: /*************************************************************************/
563: /* Check iteration limit. */
564: /*************************************************************************/
566: if (ksp->its >= max_cg_its) {
567: ksp->reason = KSP_DIVERGED_ITS;
568: PetscInfo1(ksp, "KSPSolve_STCG: iterlim: its=%d\n", ksp->its);
569: break;
570: }
572: /*************************************************************************/
573: /* Update p and the norms. */
574: /*************************************************************************/
576: VecAYPX(p, beta, z); /* p = z + beta p */
578: switch(cg->dtype) {
579: case STCG_PRECONDITIONED_DIRECTION:
580: dMp = beta*(dMp + alpha*norm_p);
581: norm_p = beta*(rzm1 + beta*norm_p);
582: break;
584: default:
585: VecDot(d, p, &dMp);
586: VecDot(p, p, &norm_p);
587: break;
588: }
590: /*************************************************************************/
591: /* Compute the new direction and update the iteration. */
592: /*************************************************************************/
594: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
595: VecDot(p, z, &kappa); /* kappa = p^T Q p */
596: ++ksp->its;
598: /*************************************************************************/
599: /* Check for negative curvature. */
600: /*************************************************************************/
602: if (kappa <= 0.0) {
603: /***********************************************************************/
604: /* In this case, the matrix is indefinite and we have encountered */
605: /* a direction of negative curvature. Follow the direction to the */
606: /* boundary of the trust region. */
607: /***********************************************************************/
609: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
610: PetscInfo1(ksp, "KSPSolve_STCG: negative curvature: kappa=%g\n", kappa);
612: if (cg->radius != 0.0 && norm_p > 0.0) {
613: /*********************************************************************/
614: /* Follow direction of negative curvature to boundary. */
615: /*********************************************************************/
617: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
618: cg->norm_d = cg->radius;
620: VecAXPY(d, step, p); /* d = d + step p */
622: /*********************************************************************/
623: /* Update objective function. */
624: /*********************************************************************/
626: cg->o_fcn += step * (0.5 * step * kappa - rz);
627: }
628: else if (cg->radius != 0.0) {
629: /*********************************************************************/
630: /* The norm of the direction is zero; there is nothing to follow. */
631: /*********************************************************************/
632: }
633: break;
634: }
635: }
637: return(0);
638: #endif
639: }
643: PetscErrorCode KSPSetUp_STCG(KSP ksp)
644: {
648: /***************************************************************************/
649: /* This implementation of CG only handles left preconditioning so generate */
650: /* an error otherwise. */
651: /***************************************************************************/
653: if (ksp->pc_side == PC_RIGHT) {
654: SETERRQ(PETSC_ERR_SUP, "No right preconditioning for KSPSTCG");
655: } else if (ksp->pc_side == PC_SYMMETRIC) {
656: SETERRQ(PETSC_ERR_SUP, "No symmetric preconditioning for KSPSTCG");
657: }
659: /***************************************************************************/
660: /* Set work vectors needed by conjugate gradient method and allocate */
661: /***************************************************************************/
663: KSPDefaultGetWork(ksp, 3);
664: return(0);
665: }
669: PetscErrorCode KSPDestroy_STCG(KSP ksp)
670: {
675: /***************************************************************************/
676: /* Clear composed functions */
677: /***************************************************************************/
679: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPSTCGSetRadius_C","",PETSC_NULL);
680: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPSTCGGetNormD_C","",PETSC_NULL);
681: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPSTCGGetObjFcn_C","",PETSC_NULL);
683: /***************************************************************************/
684: /* Destroy KSP object. */
685: /***************************************************************************/
687: KSPDefaultDestroy(ksp);
688: return(0);
689: }
694: PetscErrorCode KSPSTCGSetRadius_STCG(KSP ksp, PetscReal radius)
695: {
696: KSP_STCG *cg = (KSP_STCG *)ksp->data;
699: cg->radius = radius;
700: return(0);
701: }
705: PetscErrorCode KSPSTCGGetNormD_STCG(KSP ksp, PetscReal *norm_d)
706: {
707: KSP_STCG *cg = (KSP_STCG *)ksp->data;
710: *norm_d = cg->norm_d;
711: return(0);
712: }
716: PetscErrorCode KSPSTCGGetObjFcn_STCG(KSP ksp, PetscReal *o_fcn){
717: KSP_STCG *cg = (KSP_STCG *)ksp->data;
720: *o_fcn = cg->o_fcn;
721: return(0);
722: }
727: PetscErrorCode KSPSetFromOptions_STCG(KSP ksp)
728: {
730: KSP_STCG *cg = (KSP_STCG *)ksp->data;
733: PetscOptionsHead("KSP STCG options");
734: PetscOptionsReal("-ksp_stcg_radius", "Trust Region Radius", "KSPSTCGSetRadius", cg->radius, &cg->radius, PETSC_NULL);
735: PetscOptionsEList("-ksp_stcg_dtype", "Norm used for direction", "", DType_Table, STCG_DIRECTION_TYPES, DType_Table[cg->dtype], &cg->dtype, PETSC_NULL);
736: PetscOptionsTail();
737: return(0);
738: }
740: /*MC
741: KSPSTCG - Code to run conjugate gradient method subject to a constraint
742: on the solution norm. This is used in Trust Region methods for
743: nonlinear equations, SNESTR
745: Options Database Keys:
746: . -ksp_stcg_radius <r> - Trust Region Radius
748: Notes: This is rarely used directly
750: Use preconditioned conjugate gradient to compute
751: an approximate minimizer of the quadratic function
753: q(s) = g^T * s + 0.5 * s^T * H * s
755: subject to the trust region constraint
757: || s || <= delta,
759: where
761: delta is the trust region radius,
762: g is the gradient vector,
763: H is the Hessian approximation, and
764: M is the positive definite preconditioner matrix.
766: KSPConvergedReason may be
767: $ KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
768: $ KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
769: $ other KSP converged/diverged reasons
771: Notes:
772: The preconditioner supplied should be symmetric and positive definite.
774: Level: developer
776: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPSTCGSetRadius(), KSPSTCGGetNormD(), KSPSTCGGetObjFcn()
777: M*/
782: PetscErrorCode KSPCreate_STCG(KSP ksp)
783: {
785: KSP_STCG *cg;
789: PetscNewLog(ksp,KSP_STCG, &cg);
791: cg->radius = 0.0;
792: cg->dtype = STCG_UNPRECONDITIONED_DIRECTION;
794: ksp->data = (void *) cg;
795: ksp->pc_side = PC_LEFT;
796: ksp->normtype = KSP_NORM_UNPRECONDITIONED;
798: /***************************************************************************/
799: /* Sets the functions that are associated with this data structure */
800: /* (in C++ this is the same as defining virtual functions). */
801: /***************************************************************************/
803: ksp->ops->setup = KSPSetUp_STCG;
804: ksp->ops->solve = KSPSolve_STCG;
805: ksp->ops->destroy = KSPDestroy_STCG;
806: ksp->ops->setfromoptions = KSPSetFromOptions_STCG;
807: ksp->ops->buildsolution = KSPDefaultBuildSolution;
808: ksp->ops->buildresidual = KSPDefaultBuildResidual;
809: ksp->ops->view = 0;
811: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
812: "KSPSTCGSetRadius_C",
813: "KSPSTCGSetRadius_STCG",
814: KSPSTCGSetRadius_STCG);
815: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
816: "KSPSTCGGetNormD_C",
817: "KSPSTCGGetNormD_STCG",
818: KSPSTCGGetNormD_STCG);
819: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
820: "KSPSTCGGetObjFcn_C",
821: "KSPSTCGGetObjFcn_STCG",
822: KSPSTCGGetObjFcn_STCG);
823: return(0);
824: }