Actual source code: cross.c
1: /*
3: SLEPc singular value solver: "cross"
5: Method: Uses a Hermitian eigensolver for A^T*A
7: Last update: Jun 2007
9: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10: SLEPc - Scalable Library for Eigenvalue Problem Computations
11: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
13: This file is part of SLEPc.
14:
15: SLEPc is free software: you can redistribute it and/or modify it under the
16: terms of version 3 of the GNU Lesser General Public License as published by
17: the Free Software Foundation.
19: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
20: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
21: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
22: more details.
24: You should have received a copy of the GNU Lesser General Public License
25: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
26: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27: */
29: #include <private/svdimpl.h> /*I "slepcsvd.h" I*/
30: #include <private/epsimpl.h> /*I "slepceps.h" I*/
32: typedef struct {
33: EPS eps;
34: PetscBool setfromoptionscalled;
35: Mat mat;
36: Vec w,diag;
37: } SVD_CROSS;
41: PetscErrorCode ShellMatMult_Cross(Mat B,Vec x,Vec y)
42: {
44: SVD svd;
45: SVD_CROSS *cross;
46:
48: MatShellGetContext(B,(void**)&svd);
49: cross = (SVD_CROSS*)svd->data;
50: SVDMatMult(svd,PETSC_FALSE,x,cross->w);
51: SVDMatMult(svd,PETSC_TRUE,cross->w,y);
52: return(0);
53: }
57: PetscErrorCode ShellMatGetDiagonal_Cross(Mat B,Vec d)
58: {
59: PetscErrorCode ierr;
60: SVD svd;
61: SVD_CROSS *cross;
62: PetscInt N,n,i,j,start,end,ncols;
63: PetscScalar *work1,*work2,*diag;
64: const PetscInt *cols;
65: const PetscScalar *vals;
66:
68: MatShellGetContext(B,(void**)&svd);
69: cross = (SVD_CROSS*)svd->data;
70: if (!cross->diag) {
71: /* compute diagonal from rows and store in cross->diag */
72: VecDuplicate(d,&cross->diag);
73: SVDMatGetSize(svd,PETSC_NULL,&N);
74: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
75: PetscMalloc(sizeof(PetscScalar)*N,&work1);
76: PetscMalloc(sizeof(PetscScalar)*N,&work2);
77: for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
78: if (svd->AT) {
79: MatGetOwnershipRange(svd->AT,&start,&end);
80: for (i=start;i<end;i++) {
81: MatGetRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
82: for (j=0;j<ncols;j++)
83: work1[i] += vals[j]*vals[j];
84: MatRestoreRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
85: }
86: } else {
87: MatGetOwnershipRange(svd->A,&start,&end);
88: for (i=start;i<end;i++) {
89: MatGetRow(svd->A,i,&ncols,&cols,&vals);
90: for (j=0;j<ncols;j++)
91: work1[cols[j]] += vals[j]*vals[j];
92: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
93: }
94: }
95: MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPI_SUM,((PetscObject)svd)->comm);
96: VecGetOwnershipRange(cross->diag,&start,&end);
97: VecGetArray(cross->diag,&diag);
98: for (i=start;i<end;i++)
99: diag[i-start] = work2[i];
100: VecRestoreArray(cross->diag,&diag);
101: PetscFree(work1);
102: PetscFree(work2);
103: }
104: VecCopy(cross->diag,d);
105: return(0);
106: }
110: PetscErrorCode SVDSetUp_Cross(SVD svd)
111: {
113: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
114: PetscInt n,i;
115: PetscBool trackall;
118: if (!cross->mat) {
119: SVDMatGetLocalSize(svd,PETSC_NULL,&n);
120: MatCreateShell(((PetscObject)svd)->comm,n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
121: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))ShellMatMult_Cross);
122: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_Cross);
123: SVDMatGetVecs(svd,PETSC_NULL,&cross->w);
124: }
126: EPSSetOperators(cross->eps,cross->mat,PETSC_NULL);
127: EPSSetProblemType(cross->eps,EPS_HEP);
128: EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
129: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd);
130: EPSSetTolerances(cross->eps,svd->tol,svd->max_it);
131: /* Transfer the trackall option from svd to eps */
132: SVDGetTrackAll(svd,&trackall);
133: EPSSetTrackAll(cross->eps,trackall);
134: if (cross->setfromoptionscalled) {
135: EPSSetFromOptions(cross->eps);
136: cross->setfromoptionscalled = PETSC_FALSE;
137: }
138: EPSSetUp(cross->eps);
139: EPSGetDimensions(cross->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
140: EPSGetTolerances(cross->eps,&svd->tol,&svd->max_it);
141: /* Transfer the initial space from svd to eps */
142: if (svd->nini < 0) {
143: EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
144: for(i=0; i<-svd->nini; i++) {
145: VecDestroy(&svd->IS[i]);
146: }
147: PetscFree(svd->IS);
148: svd->nini = 0;
149: }
150: return(0);
151: }
155: PetscErrorCode SVDSolve_Cross(SVD svd)
156: {
158: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
159: PetscInt i;
160: PetscScalar sigma;
161:
163: EPSSolve(cross->eps);
164: EPSGetConverged(cross->eps,&svd->nconv);
165: EPSGetIterationNumber(cross->eps,&svd->its);
166: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
167: for (i=0;i<svd->nconv;i++) {
168: EPSGetEigenpair(cross->eps,i,&sigma,PETSC_NULL,svd->V[i],PETSC_NULL);
169: if (PetscRealPart(sigma)<0.0) SETERRQ(((PetscObject)svd)->comm,1,"Negative eigenvalue computed by EPS");
170: svd->sigma[i] = PetscSqrtReal(PetscRealPart(sigma));
171: }
172: return(0);
173: }
177: static PetscErrorCode SVDMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
178: {
179: PetscInt i;
180: SVD svd = (SVD)ctx;
181: PetscScalar er,ei;
185: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
186: er = eigr[i]; ei = eigi[i];
187: STBackTransform(eps->OP,1,&er,&ei);
188: svd->sigma[i] = PetscSqrtReal(PetscRealPart(er));
189: svd->errest[i] = errest[i];
190: }
191: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
192: return(0);
193: }
197: PetscErrorCode SVDSetFromOptions_Cross(SVD svd)
198: {
199: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
202: cross->setfromoptionscalled = PETSC_TRUE;
203: return(0);
204: }
206: EXTERN_C_BEGIN
209: PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
210: {
212: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
215: PetscObjectReference((PetscObject)eps);
216: EPSDestroy(&cross->eps);
217: cross->eps = eps;
218: PetscLogObjectParent(svd,cross->eps);
219: svd->setupcalled = 0;
220: return(0);
221: }
222: EXTERN_C_END
226: /*@
227: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
228: singular value solver.
230: Collective on SVD
232: Input Parameters:
233: + svd - singular value solver
234: - eps - the eigensolver object
236: Level: advanced
238: .seealso: SVDCrossGetEPS()
239: @*/
240: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
241: {
248: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
249: return(0);
250: }
252: EXTERN_C_BEGIN
255: PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
256: {
257: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
260: *eps = cross->eps;
261: return(0);
262: }
263: EXTERN_C_END
267: /*@
268: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
269: to the singular value solver.
271: Not Collective
273: Input Parameter:
274: . svd - singular value solver
276: Output Parameter:
277: . eps - the eigensolver object
279: Level: advanced
281: .seealso: SVDCrossSetEPS()
282: @*/
283: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
284: {
290: PetscTryMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
291: return(0);
292: }
296: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
297: {
299: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
302: PetscViewerASCIIPushTab(viewer);
303: EPSView(cross->eps,viewer);
304: PetscViewerASCIIPopTab(viewer);
305: return(0);
306: }
310: PetscErrorCode SVDReset_Cross(SVD svd)
311: {
313: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
316: EPSReset(cross->eps);
317: MatDestroy(&cross->mat);
318: VecDestroy(&cross->w);
319: VecDestroy(&cross->diag);
320: return(0);
321: }
325: PetscErrorCode SVDDestroy_Cross(SVD svd)
326: {
328: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
331: EPSDestroy(&cross->eps);
332: PetscFree(svd->data);
333: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","",PETSC_NULL);
334: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","",PETSC_NULL);
335: return(0);
336: }
338: EXTERN_C_BEGIN
341: PetscErrorCode SVDCreate_Cross(SVD svd)
342: {
344: SVD_CROSS *cross;
345: ST st;
346:
348: PetscNewLog(svd,SVD_CROSS,&cross);
349: svd->data = (void*)cross;
350: svd->ops->solve = SVDSolve_Cross;
351: svd->ops->setup = SVDSetUp_Cross;
352: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
353: svd->ops->destroy = SVDDestroy_Cross;
354: svd->ops->reset = SVDReset_Cross;
355: svd->ops->view = SVDView_Cross;
356: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","SVDCrossSetEPS_Cross",SVDCrossSetEPS_Cross);
357: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","SVDCrossGetEPS_Cross",SVDCrossGetEPS_Cross);
359: EPSCreate(((PetscObject)svd)->comm,&cross->eps);
360: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
361: EPSAppendOptionsPrefix(cross->eps,"svd_");
362: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
363: PetscLogObjectParent(svd,cross->eps);
364: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
365: EPSSetIP(cross->eps,svd->ip);
366: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
367: EPSMonitorSet(cross->eps,SVDMonitor_Cross,svd,PETSC_NULL);
368: EPSGetST(cross->eps,&st);
369: STSetMatMode(st,ST_MATMODE_SHELL);
370: cross->mat = PETSC_NULL;
371: cross->w = PETSC_NULL;
372: cross->diag = PETSC_NULL;
373: cross->setfromoptionscalled = PETSC_FALSE;
374: return(0);
375: }
376: EXTERN_C_END