xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1 
2 #include <petsc-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;
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   {
65     PetscBLASInt lierr;
66     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
67     LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr);
68     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
69     ierr = PetscFPTrapPop();CHKERRQ(ierr);
70   }
71 #else
72   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
73 #endif
74   ierr  = MatRestoreArray(jac->A,&a);CHKERRQ(ierr);
75   ierr  = MatRestoreArray(jac->U,&u);CHKERRQ(ierr);
76   ierr  = MatRestoreArray(jac->V,&v);CHKERRQ(ierr);
77   for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
78   jac->nzero = n-1-i;
79   if (jac->monitor) {
80     ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
81     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);
82     if (n >= 10) {              /* print 5 smallest and 5 largest */
83       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);
84       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);
85     } else {                    /* print all singular values */
86       char buf[256],*p;
87       size_t left = sizeof buf,used;
88       PetscInt thisline;
89       for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) {
90         ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr);
91         left -= used;
92         p += used;
93         if (thisline > 4 || i==0) {
94           ierr = PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:%s\n",buf);CHKERRQ(ierr);
95           p = buf;
96           thisline = 0;
97         }
98       }
99     }
100     ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
101   }
102   ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));
103   for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
104   for (; i<n; i++) d[i] = 0.0;
105   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);
106   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
107 #if defined(foo)
108 {
109   PetscViewer viewer;
110   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
111   ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
112   ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
113   ierr = MatView(jac->V,viewer);CHKERRQ(ierr);
114   ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
115   ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
116  }
117 #endif
118   ierr = PetscFree(work);
119   PetscFunctionReturn(0);
120 #endif
121 }
122 
123 /* -------------------------------------------------------------------------- */
124 /*
125    PCApply_SVD - Applies the SVD preconditioner to a vector.
126 
127    Input Parameters:
128 .  pc - the preconditioner context
129 .  x - input vector
130 
131    Output Parameter:
132 .  y - output vector
133 
134    Application Interface Routine: PCApply()
135  */
136 #undef __FUNCT__
137 #define __FUNCT__ "PCApply_SVD"
138 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
139 {
140   PC_SVD         *jac = (PC_SVD*)pc->data;
141   Vec            work = jac->work;
142   PetscErrorCode ierr;
143 
144   PetscFunctionBegin;
145   ierr = MatMultTranspose(jac->U,x,work);CHKERRQ(ierr);
146   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
147   ierr = MatMultTranspose(jac->V,work,y);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "PCReset_SVD"
153 static PetscErrorCode PCReset_SVD(PC pc)
154 {
155   PC_SVD         *jac = (PC_SVD*)pc->data;
156   PetscErrorCode ierr;
157 
158   PetscFunctionBegin;
159   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
160   ierr = MatDestroy(&jac->U);CHKERRQ(ierr);
161   ierr = MatDestroy(&jac->V);CHKERRQ(ierr);
162   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
163   ierr = VecDestroy(&jac->work);CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 /* -------------------------------------------------------------------------- */
168 /*
169    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
170    that was created with PCCreate_SVD().
171 
172    Input Parameter:
173 .  pc - the preconditioner context
174 
175    Application Interface Routine: PCDestroy()
176 */
177 #undef __FUNCT__
178 #define __FUNCT__ "PCDestroy_SVD"
179 static PetscErrorCode PCDestroy_SVD(PC pc)
180 {
181   PC_SVD         *jac = (PC_SVD*)pc->data;
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   ierr = PCReset_SVD(pc);CHKERRQ(ierr);
186   ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
187   ierr = PetscFree(pc->data);CHKERRQ(ierr);
188   PetscFunctionReturn(0);
189 }
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "PCSetFromOptions_SVD"
193 static PetscErrorCode PCSetFromOptions_SVD(PC pc)
194 {
195   PetscErrorCode ierr;
196   PC_SVD         *jac = (PC_SVD*)pc->data;
197   PetscBool      flg,set;
198 
199   PetscFunctionBegin;
200   ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr);
201   ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);CHKERRQ(ierr);
202   ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor?PETSC_TRUE:PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
203   if (set) {                    /* Should make PCSVDSetMonitor() */
204     if (flg && !jac->monitor) {
205       ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);CHKERRQ(ierr);
206     } else if (!flg) {
207       ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
208     }
209   }
210   ierr = PetscOptionsTail();CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 /* -------------------------------------------------------------------------- */
215 /*
216    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
217    and sets this as the private data within the generic preconditioning
218    context, PC, that was created within PCCreate().
219 
220    Input Parameter:
221 .  pc - the preconditioner context
222 
223    Application Interface Routine: PCCreate()
224 */
225 
226 /*MC
227      PCSVD - Use pseudo inverse defined by SVD of operator
228 
229    Level: advanced
230 
231   Concepts: SVD
232 
233   Options Database:
234 -  -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
235 +  -pc_svd_monitor  Print information on the extreme singular values of the operator
236 
237 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
238 M*/
239 
240 EXTERN_C_BEGIN
241 #undef __FUNCT__
242 #define __FUNCT__ "PCCreate_SVD"
243 PetscErrorCode PCCreate_SVD(PC pc)
244 {
245   PC_SVD         *jac;
246   PetscErrorCode ierr;
247 
248   PetscFunctionBegin;
249   /*
250      Creates the private data structure for this preconditioner and
251      attach it to the PC object.
252   */
253   ierr          = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr);
254   jac->zerosing = 1.e-12;
255   pc->data      = (void*)jac;
256 
257   /*
258       Set the pointers for the functions that are provided above.
259       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
260       are called, they will automatically call these functions.  Note we
261       choose not to provide a couple of these functions since they are
262       not needed.
263   */
264   pc->ops->apply               = PCApply_SVD;
265   pc->ops->applytranspose      = PCApply_SVD;
266   pc->ops->setup               = PCSetUp_SVD;
267   pc->ops->reset               = PCReset_SVD;
268   pc->ops->destroy             = PCDestroy_SVD;
269   pc->ops->setfromoptions      = PCSetFromOptions_SVD;
270   pc->ops->view                = 0;
271   pc->ops->applyrichardson     = 0;
272   PetscFunctionReturn(0);
273 }
274 EXTERN_C_END
275 
276