Actual source code: cross.c
slepc-3.5.2 2014-10-10
1: /*
3: SLEPc singular value solver: "cross"
5: Method: Uses a Hermitian eigensolver for A^T*A
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
27: #include <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
28: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
30: typedef struct {
31: EPS eps;
32: Mat mat;
33: Vec w,diag;
34: } SVD_CROSS;
38: static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
39: {
41: SVD svd;
42: SVD_CROSS *cross;
45: MatShellGetContext(B,(void**)&svd);
46: cross = (SVD_CROSS*)svd->data;
47: SVDMatMult(svd,PETSC_FALSE,x,cross->w);
48: SVDMatMult(svd,PETSC_TRUE,cross->w,y);
49: return(0);
50: }
54: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
55: {
56: PetscErrorCode ierr;
57: SVD svd;
58: SVD_CROSS *cross;
59: PetscInt N,n,i,j,start,end,ncols;
60: PetscScalar *work1,*work2,*diag;
61: const PetscInt *cols;
62: const PetscScalar *vals;
65: MatShellGetContext(B,(void**)&svd);
66: cross = (SVD_CROSS*)svd->data;
67: if (!cross->diag) {
68: /* compute diagonal from rows and store in cross->diag */
69: VecDuplicate(d,&cross->diag);
70: SVDMatGetSize(svd,NULL,&N);
71: SVDMatGetLocalSize(svd,NULL,&n);
72: PetscMalloc2(N,&work1,N,&work2);
73: for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
74: if (svd->AT) {
75: MatGetOwnershipRange(svd->AT,&start,&end);
76: for (i=start;i<end;i++) {
77: MatGetRow(svd->AT,i,&ncols,NULL,&vals);
78: for (j=0;j<ncols;j++)
79: work1[i] += vals[j]*vals[j];
80: MatRestoreRow(svd->AT,i,&ncols,NULL,&vals);
81: }
82: } else {
83: MatGetOwnershipRange(svd->A,&start,&end);
84: for (i=start;i<end;i++) {
85: MatGetRow(svd->A,i,&ncols,&cols,&vals);
86: for (j=0;j<ncols;j++)
87: work1[cols[j]] += vals[j]*vals[j];
88: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
89: }
90: }
91: MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)svd));
92: VecGetOwnershipRange(cross->diag,&start,&end);
93: VecGetArray(cross->diag,&diag);
94: for (i=start;i<end;i++)
95: diag[i-start] = work2[i];
96: VecRestoreArray(cross->diag,&diag);
97: PetscFree2(work1,work2);
98: }
99: VecCopy(cross->diag,d);
100: return(0);
101: }
105: PetscErrorCode SVDSetUp_Cross(SVD svd)
106: {
108: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
109: PetscInt n;
110: PetscBool trackall;
113: if (!cross->mat) {
114: SVDMatGetLocalSize(svd,NULL,&n);
115: MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
116: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
117: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
118: SVDMatGetVecs(svd,NULL,&cross->w);
119: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->mat);
120: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->w);
121: }
123: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
124: EPSSetOperators(cross->eps,cross->mat,NULL);
125: EPSSetProblemType(cross->eps,EPS_HEP);
126: EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
127: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv?svd->ncv:PETSC_DEFAULT,svd->mpd?svd->mpd:PETSC_DEFAULT);
128: EPSSetTolerances(cross->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it?svd->max_it:PETSC_DEFAULT);
129: /* Transfer the trackall option from svd to eps */
130: SVDGetTrackAll(svd,&trackall);
131: EPSSetTrackAll(cross->eps,trackall);
132: EPSSetUp(cross->eps);
133: EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
134: EPSGetTolerances(cross->eps,NULL,&svd->max_it);
135: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
136: /* Transfer the initial space from svd to eps */
137: if (svd->nini < 0) {
138: EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
139: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
140: }
141: svd->leftbasis = PETSC_FALSE;
142: SVDAllocateSolution(svd,0);
143: return(0);
144: }
148: PetscErrorCode SVDSolve_Cross(SVD svd)
149: {
151: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
152: PetscInt i;
153: PetscScalar sigma;
154: Vec v;
157: EPSSolve(cross->eps);
158: EPSGetConverged(cross->eps,&svd->nconv);
159: EPSGetIterationNumber(cross->eps,&svd->its);
160: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
161: for (i=0;i<svd->nconv;i++) {
162: BVGetColumn(svd->V,i,&v);
163: EPSGetEigenpair(cross->eps,i,&sigma,NULL,v,NULL);
164: BVRestoreColumn(svd->V,i,&v);
165: if (PetscRealPart(sigma)<0.0) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS");
166: svd->sigma[i] = PetscSqrtReal(PetscRealPart(sigma));
167: }
168: return(0);
169: }
173: static PetscErrorCode SVDMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
174: {
175: PetscInt i;
176: SVD svd = (SVD)ctx;
177: PetscScalar er,ei;
181: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
182: er = eigr[i]; ei = eigi[i];
183: STBackTransform(eps->st,1,&er,&ei);
184: svd->sigma[i] = PetscSqrtReal(PetscRealPart(er));
185: svd->errest[i] = errest[i];
186: }
187: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
188: return(0);
189: }
193: PetscErrorCode SVDSetFromOptions_Cross(SVD svd)
194: {
196: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
199: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
200: EPSSetFromOptions(cross->eps);
201: return(0);
202: }
206: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
207: {
209: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
212: PetscObjectReference((PetscObject)eps);
213: EPSDestroy(&cross->eps);
214: cross->eps = eps;
215: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
216: svd->setupcalled = 0;
217: return(0);
218: }
222: /*@
223: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
224: singular value solver.
226: Collective on SVD
228: Input Parameters:
229: + svd - singular value solver
230: - eps - the eigensolver object
232: Level: advanced
234: .seealso: SVDCrossGetEPS()
235: @*/
236: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
237: {
244: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
245: return(0);
246: }
250: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
251: {
252: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
253: ST st;
257: if (!cross->eps) {
258: EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
259: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
260: EPSAppendOptionsPrefix(cross->eps,"svd_");
261: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
262: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
263: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
264: EPSMonitorSet(cross->eps,SVDMonitor_Cross,svd,NULL);
265: EPSGetST(cross->eps,&st);
266: STSetMatMode(st,ST_MATMODE_SHELL);
267: }
268: *eps = cross->eps;
269: return(0);
270: }
274: /*@
275: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
276: to the singular value solver.
278: Not Collective
280: Input Parameter:
281: . svd - singular value solver
283: Output Parameter:
284: . eps - the eigensolver object
286: Level: advanced
288: .seealso: SVDCrossSetEPS()
289: @*/
290: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
291: {
297: PetscTryMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
298: return(0);
299: }
303: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
304: {
306: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
309: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
310: PetscViewerASCIIPushTab(viewer);
311: EPSView(cross->eps,viewer);
312: PetscViewerASCIIPopTab(viewer);
313: return(0);
314: }
318: PetscErrorCode SVDReset_Cross(SVD svd)
319: {
321: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
324: if (cross->eps) { EPSReset(cross->eps); }
325: MatDestroy(&cross->mat);
326: VecDestroy(&cross->w);
327: VecDestroy(&cross->diag);
328: return(0);
329: }
333: PetscErrorCode SVDDestroy_Cross(SVD svd)
334: {
336: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
339: EPSDestroy(&cross->eps);
340: PetscFree(svd->data);
341: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
342: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
343: return(0);
344: }
348: PETSC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
349: {
351: SVD_CROSS *cross;
354: PetscNewLog(svd,&cross);
355: svd->data = (void*)cross;
357: svd->ops->solve = SVDSolve_Cross;
358: svd->ops->setup = SVDSetUp_Cross;
359: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
360: svd->ops->destroy = SVDDestroy_Cross;
361: svd->ops->reset = SVDReset_Cross;
362: svd->ops->view = SVDView_Cross;
363: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
364: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
365: return(0);
366: }