xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 7cf18bfc0af3bd20f57ca1ab9dbfc3d774cdee5d)
1 
2 #include <private/pcimpl.h>   /*I "petscpc.h" I*/
3 #include <petscblaslapack.h>
4 
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;
15 
16 
17 /* -------------------------------------------------------------------------- */
18 /*
19    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
20                     by setting data structures and options.
21 
22    Input Parameter:
23 .  pc - the preconditioner context
24 
25    Application Interface Routine: PCSetUp()
26 
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 */
31 #undef __FUNCT__
32 #define __FUNCT__ "PCSetUp_SVD"
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;
39   PetscErrorCode ierr;
40   PetscScalar    *a,*u,*v,*d,*work;
41   PetscBLASInt   nb,lwork,lierr;
42   PetscInt       i,n;
43 
44   PetscFunctionBegin;
45   if (!jac->diag) {
46     /* assume square matrices */
47     ierr = MatGetVecs(pc->pmat,&jac->diag,&jac->work);CHKERRQ(ierr);
48   }
49   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
50   ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
51   if (!jac->U) {
52     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr);
53     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->V);CHKERRQ(ierr);
54   }
55   ierr = MatGetSize(pc->pmat,&n,PETSC_NULL);CHKERRQ(ierr);
56   nb    = PetscBLASIntCast(n);
57   lwork = 5*nb;
58   ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
59   ierr  = MatGetArray(jac->A,&a);CHKERRQ(ierr);
60   ierr  = MatGetArray(jac->U,&u);CHKERRQ(ierr);
61   ierr  = MatGetArray(jac->V,&v);CHKERRQ(ierr);
62   ierr  = VecGetArray(jac->diag,&d);CHKERRQ(ierr);
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   ierr  = MatRestoreArray(jac->A,&a);CHKERRQ(ierr);
70   ierr  = MatRestoreArray(jac->U,&u);CHKERRQ(ierr);
71   ierr  = MatRestoreArray(jac->V,&v);CHKERRQ(ierr);
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     ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
76     ierr = 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);CHKERRQ(ierr);
77     if (n >= 10) {              /* print 5 smallest and 5 largest */
78       ierr = 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]));CHKERRQ(ierr);
79       ierr = 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]));CHKERRQ(ierr);
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         ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr);
86         left -= used;
87         p += used;
88         if (thisline > 4 || i==0) {
89           ierr = PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:%s\n",buf);CHKERRQ(ierr);
90           p = buf;
91           thisline = 0;
92         }
93       }
94     }
95     ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
96   }
97   ierr = 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   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);
101   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
102 #if defined(foo)
103 {
104   PetscViewer viewer;
105   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
106   ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
107   ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
108   ierr = MatView(jac->V,viewer);CHKERRQ(ierr);
109   ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
110   ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
111  }
112 #endif
113   ierr = PetscFree(work);
114   PetscFunctionReturn(0);
115 #endif
116 }
117 
118 /* -------------------------------------------------------------------------- */
119 /*
120    PCApply_SVD - Applies the SVD preconditioner to a vector.
121 
122    Input Parameters:
123 .  pc - the preconditioner context
124 .  x - input vector
125 
126    Output Parameter:
127 .  y - output vector
128 
129    Application Interface Routine: PCApply()
130  */
131 #undef __FUNCT__
132 #define __FUNCT__ "PCApply_SVD"
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;
137   PetscErrorCode ierr;
138 
139   PetscFunctionBegin;
140   ierr = MatMultTranspose(jac->U,x,work);CHKERRQ(ierr);
141   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
142   ierr = MatMultTranspose(jac->V,work,y);CHKERRQ(ierr);
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "PCReset_SVD"
148 static PetscErrorCode PCReset_SVD(PC pc)
149 {
150   PC_SVD         *jac = (PC_SVD*)pc->data;
151   PetscErrorCode ierr;
152 
153   PetscFunctionBegin;
154   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
155   ierr = MatDestroy(&jac->U);CHKERRQ(ierr);
156   ierr = MatDestroy(&jac->V);CHKERRQ(ierr);
157   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
158   ierr = VecDestroy(&jac->work);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 /* -------------------------------------------------------------------------- */
163 /*
164    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
165    that was created with PCCreate_SVD().
166 
167    Input Parameter:
168 .  pc - the preconditioner context
169 
170    Application Interface Routine: PCDestroy()
171 */
172 #undef __FUNCT__
173 #define __FUNCT__ "PCDestroy_SVD"
174 static PetscErrorCode PCDestroy_SVD(PC pc)
175 {
176   PC_SVD         *jac = (PC_SVD*)pc->data;
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   ierr = PCReset_SVD(pc);CHKERRQ(ierr);
181   ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
182   ierr = PetscFree(pc->data);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PCSetFromOptions_SVD"
188 static PetscErrorCode PCSetFromOptions_SVD(PC pc)
189 {
190   PetscErrorCode ierr;
191   PC_SVD         *jac = (PC_SVD*)pc->data;
192   PetscBool      flg,set;
193 
194   PetscFunctionBegin;
195   ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr);
196   ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);CHKERRQ(ierr);
197   ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremeal singular values","None",jac->monitor?PETSC_TRUE:PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
198   if (set) {                    /* Should make PCSVDSetMonitor() */
199     if (flg && !jac->monitor) {
200       ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);CHKERRQ(ierr);
201     } else if (!flg) {
202       ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
203     }
204   }
205   ierr = PetscOptionsTail();CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
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().
214 
215    Input Parameter:
216 .  pc - the preconditioner context
217 
218    Application Interface Routine: PCCreate()
219 */
220 
221 /*MC
222      PCSVD - Use pseudo inverse defined by SVD of operator
223 
224    Level: advanced
225 
226   Concepts: SVD
227 
228   Options Database:
229 .  -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
230 
231 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
232 M*/
233 
234 EXTERN_C_BEGIN
235 #undef __FUNCT__
236 #define __FUNCT__ "PCCreate_SVD"
237 PetscErrorCode PCCreate_SVD(PC pc)
238 {
239   PC_SVD         *jac;
240   PetscErrorCode ierr;
241 
242   PetscFunctionBegin;
243   /*
244      Creates the private data structure for this preconditioner and
245      attach it to the PC object.
246   */
247   ierr          = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr);
248   jac->zerosing = 1.e-12;
249   pc->data      = (void*)jac;
250 
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   PetscFunctionReturn(0);
267 }
268 EXTERN_C_END
269 
270