Actual source code: cross.c
slepc-3.6.1 2015-09-03
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-2015, 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++) diag[i-start] = work2[i];
95: VecRestoreArray(cross->diag,&diag);
96: PetscFree2(work1,work2);
97: }
98: VecCopy(cross->diag,d);
99: return(0);
100: }
104: PetscErrorCode SVDSetUp_Cross(SVD svd)
105: {
107: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
108: PetscInt n;
109: PetscBool trackall;
112: if (!cross->mat) {
113: SVDMatGetLocalSize(svd,NULL,&n);
114: MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
115: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
116: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
117: SVDMatCreateVecs(svd,NULL,&cross->w);
118: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->mat);
119: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->w);
120: }
122: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
123: EPSSetOperators(cross->eps,cross->mat,NULL);
124: EPSSetProblemType(cross->eps,EPS_HEP);
125: EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
126: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv?svd->ncv:PETSC_DEFAULT,svd->mpd?svd->mpd:PETSC_DEFAULT);
127: EPSSetTolerances(cross->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it?svd->max_it:PETSC_DEFAULT);
128: /* Transfer the trackall option from svd to eps */
129: SVDGetTrackAll(svd,&trackall);
130: EPSSetTrackAll(cross->eps,trackall);
131: EPSSetUp(cross->eps);
132: EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
133: EPSGetTolerances(cross->eps,NULL,&svd->max_it);
134: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
135: /* Transfer the initial space from svd to eps */
136: if (svd->nini < 0) {
137: EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
138: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
139: }
140: svd->leftbasis = PETSC_FALSE;
141: SVDAllocateSolution(svd,0);
142: return(0);
143: }
147: PetscErrorCode SVDSolve_Cross(SVD svd)
148: {
150: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
151: PetscInt i;
152: PetscScalar sigma;
153: Vec v;
156: EPSSolve(cross->eps);
157: EPSGetConverged(cross->eps,&svd->nconv);
158: EPSGetIterationNumber(cross->eps,&svd->its);
159: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
160: for (i=0;i<svd->nconv;i++) {
161: BVGetColumn(svd->V,i,&v);
162: EPSGetEigenpair(cross->eps,i,&sigma,NULL,v,NULL);
163: BVRestoreColumn(svd->V,i,&v);
164: if (PetscRealPart(sigma)<0.0) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS");
165: svd->sigma[i] = PetscSqrtReal(PetscRealPart(sigma));
166: }
167: return(0);
168: }
172: static PetscErrorCode SVDMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
173: {
174: PetscInt i;
175: SVD svd = (SVD)ctx;
176: PetscScalar er,ei;
180: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
181: er = eigr[i]; ei = eigi[i];
182: STBackTransform(eps->st,1,&er,&ei);
183: svd->sigma[i] = PetscSqrtReal(PetscRealPart(er));
184: svd->errest[i] = errest[i];
185: }
186: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
187: return(0);
188: }
192: PetscErrorCode SVDSetFromOptions_Cross(PetscOptions *PetscOptionsObject,SVD svd)
193: {
195: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
198: PetscOptionsHead(PetscOptionsObject,"SVD Cross Options");
199: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
200: EPSSetFromOptions(cross->eps);
201: PetscOptionsTail();
202: return(0);
203: }
207: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
208: {
210: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
213: PetscObjectReference((PetscObject)eps);
214: EPSDestroy(&cross->eps);
215: cross->eps = eps;
216: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
217: svd->state = SVD_STATE_INITIAL;
218: return(0);
219: }
223: /*@
224: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
225: singular value solver.
227: Collective on SVD
229: Input Parameters:
230: + svd - singular value solver
231: - eps - the eigensolver object
233: Level: advanced
235: .seealso: SVDCrossGetEPS()
236: @*/
237: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
238: {
245: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
246: return(0);
247: }
251: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
252: {
253: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
254: ST st;
258: if (!cross->eps) {
259: EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
260: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
261: EPSAppendOptionsPrefix(cross->eps,"svd_");
262: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
263: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
264: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
265: EPSMonitorSet(cross->eps,SVDMonitor_Cross,svd,NULL);
266: EPSGetST(cross->eps,&st);
267: STSetMatMode(st,ST_MATMODE_SHELL);
268: }
269: *eps = cross->eps;
270: return(0);
271: }
275: /*@
276: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
277: to the singular value solver.
279: Not Collective
281: Input Parameter:
282: . svd - singular value solver
284: Output Parameter:
285: . eps - the eigensolver object
287: Level: advanced
289: .seealso: SVDCrossSetEPS()
290: @*/
291: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
292: {
298: PetscTryMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
299: return(0);
300: }
304: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
305: {
307: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
310: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
311: PetscViewerASCIIPushTab(viewer);
312: EPSView(cross->eps,viewer);
313: PetscViewerASCIIPopTab(viewer);
314: return(0);
315: }
319: PetscErrorCode SVDReset_Cross(SVD svd)
320: {
322: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
325: if (cross->eps) { EPSReset(cross->eps); }
326: MatDestroy(&cross->mat);
327: VecDestroy(&cross->w);
328: VecDestroy(&cross->diag);
329: return(0);
330: }
334: PetscErrorCode SVDDestroy_Cross(SVD svd)
335: {
337: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
340: EPSDestroy(&cross->eps);
341: PetscFree(svd->data);
342: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
343: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
344: return(0);
345: }
349: PETSC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
350: {
352: SVD_CROSS *cross;
355: PetscNewLog(svd,&cross);
356: svd->data = (void*)cross;
358: svd->ops->solve = SVDSolve_Cross;
359: svd->ops->setup = SVDSetUp_Cross;
360: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
361: svd->ops->destroy = SVDDestroy_Cross;
362: svd->ops->reset = SVDReset_Cross;
363: svd->ops->view = SVDView_Cross;
364: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
365: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
366: return(0);
367: }