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