Actual source code: svd.c
2: #include <private/pcimpl.h> /*I "petscpc.h" I*/
3: #include <petscblaslapack.h>
5: /*
6: Private context (data structure) for the SVD preconditioner.
7: */
8: typedef struct {
9: Vec diag,work;
10: Mat A,U,V;
11: PetscInt nzero;
12: PetscReal zerosing; /* measure of smallest singular value treated as nonzero */
13: PetscViewer monitor;
14: } PC_SVD;
17: /* -------------------------------------------------------------------------- */
18: /*
19: PCSetUp_SVD - Prepares for the use of the SVD preconditioner
20: by setting data structures and options.
22: Input Parameter:
23: . pc - the preconditioner context
25: Application Interface Routine: PCSetUp()
27: Notes:
28: The interface routine PCSetUp() is not usually called directly by
29: the user, but instead is called by PCApply() if necessary.
30: */
33: static PetscErrorCode PCSetUp_SVD(PC pc)
34: {
35: #if defined(PETSC_MISSING_LAPACK_GESVD)
36: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates.");
37: #else
38: PC_SVD *jac = (PC_SVD*)pc->data;
40: PetscScalar *a,*u,*v,*d,*work;
41: PetscBLASInt nb,lwork,lierr;
42: PetscInt i,n;
45: if (!jac->diag) {
46: /* assume square matrices */
47: MatGetVecs(pc->pmat,&jac->diag,&jac->work);
48: }
49: MatDestroy(&jac->A);
50: MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);
51: if (!jac->U) {
52: MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);
53: MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->V);
54: }
55: MatGetSize(pc->pmat,&n,PETSC_NULL);
56: nb = PetscBLASIntCast(n);
57: lwork = 5*nb;
58: PetscMalloc(lwork*sizeof(PetscScalar),&work);
59: MatGetArray(jac->A,&a);
60: MatGetArray(jac->U,&u);
61: MatGetArray(jac->V,&v);
62: VecGetArray(jac->diag,&d);
63: #if !defined(PETSC_USE_COMPLEX)
64: LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr);
65: #else
66: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
67: #endif
68: if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
69: MatRestoreArray(jac->A,&a);
70: MatRestoreArray(jac->U,&u);
71: MatRestoreArray(jac->V,&v);
72: for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
73: jac->nzero = n-1-i;
74: if (jac->monitor) {
75: PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);
76: PetscViewerASCIIPrintf(jac->monitor," SVD: condition number %14.12e, %D of %D singular values are (nearly) zero\n",(double)PetscRealPart(d[0]/d[n-1]),jac->nzero,n);
77: if (n >= 10) { /* print 5 smallest and 5 largest */
78: PetscViewerASCIIPrintf(jac->monitor," SVD: smallest singular values: %14.12e %14.12e %14.12e %14.12e %14.12e\n",(double)PetscRealPart(d[n-1]),(double)PetscRealPart(d[n-2]),(double)PetscRealPart(d[n-3]),(double)PetscRealPart(d[n-4]),(double)PetscRealPart(d[n-5]));
79: PetscViewerASCIIPrintf(jac->monitor," SVD: largest singular values : %14.12e %14.12e %14.12e %14.12e %14.12e\n",(double)PetscRealPart(d[4]),(double)PetscRealPart(d[3]),(double)PetscRealPart(d[2]),(double)PetscRealPart(d[1]),(double)PetscRealPart(d[0]));
80: } else { /* print all singular values */
81: char buf[256],*p;
82: size_t left = sizeof buf,used;
83: PetscInt thisline;
84: for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) {
85: PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));
86: left -= used;
87: p += used;
88: if (thisline > 4 || i==0) {
89: PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);
90: p = buf;
91: thisline = 0;
92: }
93: }
94: }
95: PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);
96: }
97: PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));
98: for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
99: for (; i<n; i++) d[i] = 0.0;
100: PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);
101: VecRestoreArray(jac->diag,&d);
102: #if defined(foo)
103: {
104: PetscViewer viewer;
105: PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);
106: MatView(jac->A,viewer);
107: MatView(jac->U,viewer);
108: MatView(jac->V,viewer);
109: VecView(jac->diag,viewer);
110: PetscViewerDestroy(viewer);
111: }
112: #endif
113: PetscFree(work);
114: return(0);
115: #endif
116: }
118: /* -------------------------------------------------------------------------- */
119: /*
120: PCApply_SVD - Applies the SVD preconditioner to a vector.
122: Input Parameters:
123: . pc - the preconditioner context
124: . x - input vector
126: Output Parameter:
127: . y - output vector
129: Application Interface Routine: PCApply()
130: */
133: static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
134: {
135: PC_SVD *jac = (PC_SVD*)pc->data;
136: Vec work = jac->work;
140: MatMultTranspose(jac->U,x,work);
141: VecPointwiseMult(work,work,jac->diag);
142: MatMultTranspose(jac->V,work,y);
143: return(0);
144: }
148: static PetscErrorCode PCReset_SVD(PC pc)
149: {
150: PC_SVD *jac = (PC_SVD*)pc->data;
154: MatDestroy(&jac->A);
155: MatDestroy(&jac->U);
156: MatDestroy(&jac->V);
157: VecDestroy(&jac->diag);
158: VecDestroy(&jac->work);
159: return(0);
160: }
162: /* -------------------------------------------------------------------------- */
163: /*
164: PCDestroy_SVD - Destroys the private context for the SVD preconditioner
165: that was created with PCCreate_SVD().
167: Input Parameter:
168: . pc - the preconditioner context
170: Application Interface Routine: PCDestroy()
171: */
174: static PetscErrorCode PCDestroy_SVD(PC pc)
175: {
176: PC_SVD *jac = (PC_SVD*)pc->data;
180: PCReset_SVD(pc);
181: PetscViewerDestroy(&jac->monitor);
182: PetscFree(pc->data);
183: return(0);
184: }
188: static PetscErrorCode PCSetFromOptions_SVD(PC pc)
189: {
191: PC_SVD *jac = (PC_SVD*)pc->data;
192: PetscBool flg,set;
195: PetscOptionsHead("SVD options");
196: PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);
197: PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremeal singular values","None",jac->monitor?PETSC_TRUE:PETSC_FALSE,&flg,&set);
198: if (set) { /* Should make PCSVDSetMonitor() */
199: if (flg && !jac->monitor) {
200: PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);
201: } else if (!flg) {
202: PetscViewerDestroy(&jac->monitor);
203: }
204: }
205: PetscOptionsTail();
206: return(0);
207: }
209: /* -------------------------------------------------------------------------- */
210: /*
211: PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
212: and sets this as the private data within the generic preconditioning
213: context, PC, that was created within PCCreate().
215: Input Parameter:
216: . pc - the preconditioner context
218: Application Interface Routine: PCCreate()
219: */
221: /*MC
222: PCSVD - Use pseudo inverse defined by SVD of operator
224: Level: advanced
226: Concepts: SVD
228: Options Database:
229: . -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
231: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
232: M*/
237: PetscErrorCode PCCreate_SVD(PC pc)
238: {
239: PC_SVD *jac;
243: /*
244: Creates the private data structure for this preconditioner and
245: attach it to the PC object.
246: */
247: PetscNewLog(pc,PC_SVD,&jac);
248: jac->zerosing = 1.e-12;
249: pc->data = (void*)jac;
251: /*
252: Set the pointers for the functions that are provided above.
253: Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
254: are called, they will automatically call these functions. Note we
255: choose not to provide a couple of these functions since they are
256: not needed.
257: */
258: pc->ops->apply = PCApply_SVD;
259: pc->ops->applytranspose = PCApply_SVD;
260: pc->ops->setup = PCSetUp_SVD;
261: pc->ops->reset = PCReset_SVD;
262: pc->ops->destroy = PCDestroy_SVD;
263: pc->ops->setfromoptions = PCSetFromOptions_SVD;
264: pc->ops->view = 0;
265: pc->ops->applyrichardson = 0;
266: return(0);
267: }