xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
1 #define PETSCKSP_DLL
2 
3 #include "private/pcimpl.h"   /*I "petscpc.h" I*/
4 #include "petscblaslapack.h"
5 
6 /*
7    Private context (data structure) for the SVD preconditioner.
8 */
9 typedef struct {
10   Vec        diag,work;
11   Mat        A,U,V;
12   PetscInt   nzero;
13 } PC_SVD;
14 
15 
16 /* -------------------------------------------------------------------------- */
17 /*
18    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
19                     by setting data structures and options.
20 
21    Input Parameter:
22 .  pc - the preconditioner context
23 
24    Application Interface Routine: PCSetUp()
25 
26    Notes:
27    The interface routine PCSetUp() is not usually called directly by
28    the user, but instead is called by PCApply() if necessary.
29 */
30 #undef __FUNCT__
31 #define __FUNCT__ "PCSetUp_SVD"
32 static PetscErrorCode PCSetUp_SVD(PC pc)
33 {
34   PC_SVD         *jac = (PC_SVD*)pc->data;
35   PetscErrorCode ierr;
36   PetscScalar    *a,*u,*v,*d,*work;
37   PetscBLASInt   nb,lwork,lierr;
38   PetscInt       i,n;
39 
40   PetscFunctionBegin;
41   if (!jac->diag) {
42     /* assume square matrices */
43     ierr = MatGetVecs(pc->mat,&jac->diag,&jac->work);CHKERRQ(ierr);
44   }
45   if (jac->A) {
46     ierr = MatDestroy(jac->A);CHKERRQ(ierr);
47   }
48   ierr = MatConvert(pc->mat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
49   if (!jac->U) {
50     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr);
51     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->V);CHKERRQ(ierr);
52   }
53   ierr = MatGetSize(pc->mat,&n,PETSC_NULL);CHKERRQ(ierr);
54   nb    = PetscBLASIntCast(n);
55   lwork = 5*nb;
56   ierr  = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
57   ierr  = MatGetArray(jac->A,&a);CHKERRQ(ierr);
58   ierr  = MatGetArray(jac->U,&u);CHKERRQ(ierr);
59   ierr  = MatGetArray(jac->V,&v);CHKERRQ(ierr);
60   ierr  = VecGetArray(jac->diag,&d);CHKERRQ(ierr);
61 #if !defined(PETSC_USE_COMPLEX)
62   LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr);
63 #else
64   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
65 #endif
66   if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
67   ierr  = MatRestoreArray(jac->A,&a);CHKERRQ(ierr);
68   ierr  = MatRestoreArray(jac->U,&u);CHKERRQ(ierr);
69   ierr  = MatRestoreArray(jac->V,&v);CHKERRQ(ierr);
70   jac->nzero = 0;
71   for (i=0; i<n; i++) {
72     if (PetscRealPart(d[i]) < 1.e-12) {jac->nzero = n - i;break;}
73     d[i] = 1.0/d[i];
74   }
75   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);
76   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
77 #if defined(foo)
78 {
79   PetscViewer viewer;
80   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
81   ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
82   ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
83   ierr = MatView(jac->V,viewer);CHKERRQ(ierr);
84   ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
85   ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
86  }
87 #endif
88   ierr = PetscFree(work);
89   PetscFunctionReturn(0);
90 }
91 
92 /* -------------------------------------------------------------------------- */
93 /*
94    PCApply_SVD - Applies the SVD preconditioner to a vector.
95 
96    Input Parameters:
97 .  pc - the preconditioner context
98 .  x - input vector
99 
100    Output Parameter:
101 .  y - output vector
102 
103    Application Interface Routine: PCApply()
104  */
105 #undef __FUNCT__
106 #define __FUNCT__ "PCApply_SVD"
107 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
108 {
109   PC_SVD         *jac = (PC_SVD*)pc->data;
110   Vec            work = jac->work;
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114   ierr = MatMultTranspose(jac->U,x,work);CHKERRQ(ierr);
115   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
116   ierr = MatMultTranspose(jac->V,work,y);CHKERRQ(ierr);
117   PetscFunctionReturn(0);
118 }
119 
120 /* -------------------------------------------------------------------------- */
121 /*
122    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
123    that was created with PCCreate_SVD().
124 
125    Input Parameter:
126 .  pc - the preconditioner context
127 
128    Application Interface Routine: PCDestroy()
129 */
130 #undef __FUNCT__
131 #define __FUNCT__ "PCDestroy_SVD"
132 static PetscErrorCode PCDestroy_SVD(PC pc)
133 {
134   PC_SVD         *jac = (PC_SVD*)pc->data;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   if (jac->A) {
139     ierr = MatDestroy(jac->A);CHKERRQ(ierr);
140   }
141   if (jac->U) {
142     ierr = MatDestroy(jac->U);CHKERRQ(ierr);
143   }
144   if (jac->V) {
145     ierr = MatDestroy(jac->V);CHKERRQ(ierr);
146   }
147   if (jac->diag) {
148     ierr = VecDestroy(jac->diag);CHKERRQ(ierr);
149   }
150   ierr = PetscFree(jac);CHKERRQ(ierr);
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNCT__
155 #define __FUNCT__ "PCSetFromOptions_SVD"
156 static PetscErrorCode PCSetFromOptions_SVD(PC pc)
157 {
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr);
162   ierr = PetscOptionsTail();CHKERRQ(ierr);
163   PetscFunctionReturn(0);
164 }
165 
166 /* -------------------------------------------------------------------------- */
167 /*
168    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
169    and sets this as the private data within the generic preconditioning
170    context, PC, that was created within PCCreate().
171 
172    Input Parameter:
173 .  pc - the preconditioner context
174 
175    Application Interface Routine: PCCreate()
176 */
177 
178 /*MC
179      PCSVD - Use pseudo inverse defined by SVD of operator
180 
181    Level: advanced
182 
183   Concepts: SVD
184 
185          Zero entries along the diagonal are replaced with the value 0.0
186 
187 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
188 M*/
189 
190 EXTERN_C_BEGIN
191 #undef __FUNCT__
192 #define __FUNCT__ "PCCreate_SVD"
193 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SVD(PC pc)
194 {
195   PC_SVD         *jac;
196   PetscErrorCode ierr;
197 
198   PetscFunctionBegin;
199   /*
200      Creates the private data structure for this preconditioner and
201      attach it to the PC object.
202   */
203   ierr      = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr);
204   pc->data  = (void*)jac;
205 
206   /*
207       Set the pointers for the functions that are provided above.
208       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
209       are called, they will automatically call these functions.  Note we
210       choose not to provide a couple of these functions since they are
211       not needed.
212   */
213   pc->ops->apply               = PCApply_SVD;
214   pc->ops->applytranspose      = PCApply_SVD;
215   pc->ops->setup               = PCSetUp_SVD;
216   pc->ops->destroy             = PCDestroy_SVD;
217   pc->ops->setfromoptions      = PCSetFromOptions_SVD;
218   pc->ops->view                = 0;
219   pc->ops->applyrichardson     = 0;
220   PetscFunctionReturn(0);
221 }
222 EXTERN_C_END
223 
224