Actual source code: cyclic.c
1: /*
3: SLEPc singular value solver: "cyclic"
5: Method: Uses a Hermitian eigensolver for H(A) = [ 0 A ; A^T 0 ]
7: Last update: Jun 2007
9: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10: SLEPc - Scalable Library for Eigenvalue Problem Computations
11: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
13: This file is part of SLEPc.
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 <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
30: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
32: typedef struct {
33: PetscBool explicitmatrix;
34: EPS eps;
35: PetscBool setfromoptionscalled;
36: Mat mat;
37: Vec x1,x2,y1,y2;
38: } SVD_CYCLIC;
42: static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
43: {
44: PetscErrorCode ierr;
45: SVD svd;
46: SVD_CYCLIC *cyclic;
47: const PetscScalar *px;
48: PetscScalar *py;
49: PetscInt m;
52: MatShellGetContext(B,(void**)&svd);
53: cyclic = (SVD_CYCLIC*)svd->data;
54: SVDMatGetLocalSize(svd,&m,NULL);
55: VecGetArrayRead(x,&px);
56: VecGetArray(y,&py);
57: VecPlaceArray(cyclic->x1,px);
58: VecPlaceArray(cyclic->x2,px+m);
59: VecPlaceArray(cyclic->y1,py);
60: VecPlaceArray(cyclic->y2,py+m);
61: SVDMatMult(svd,PETSC_FALSE,cyclic->x2,cyclic->y1);
62: SVDMatMult(svd,PETSC_TRUE,cyclic->x1,cyclic->y2);
63: VecResetArray(cyclic->x1);
64: VecResetArray(cyclic->x2);
65: VecResetArray(cyclic->y1);
66: VecResetArray(cyclic->y2);
67: VecRestoreArrayRead(x,&px);
68: VecRestoreArray(y,&py);
69: return(0);
70: }
74: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
75: {
79: VecSet(diag,0.0);
80: return(0);
81: }
85: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
86: {
87: PetscErrorCode ierr;
88: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
89: PetscInt M,N,m,n,i,isl,Istart,Iend;
90: const PetscScalar *isa;
91: PetscScalar *va;
92: PetscBool trackall;
93: Vec v;
94: Mat Zm,Zn;
97: SVDMatGetSize(svd,&M,&N);
98: SVDMatGetLocalSize(svd,&m,&n);
99: if (!cyclic->mat) {
100: if (cyclic->explicitmatrix) {
101: if (!svd->AT) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
102: MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
103: MatSetSizes(Zm,m,m,M,M);
104: MatSetFromOptions(Zm);
105: MatSetUp(Zm);
106: MatGetOwnershipRange(Zm,&Istart,&Iend);
107: for (i=Istart;i<Iend;i++) {
108: MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
109: }
110: MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
111: MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
112: MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
113: MatSetSizes(Zn,n,n,N,N);
114: MatSetFromOptions(Zn);
115: MatSetUp(Zn);
116: MatGetOwnershipRange(Zn,&Istart,&Iend);
117: for (i=Istart;i<Iend;i++) {
118: MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
119: }
120: MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
121: MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
122: SlepcMatTile(1.0,Zm,1.0,svd->A,1.0,svd->AT,1.0,Zn,&cyclic->mat);
123: PetscLogObjectParent(svd,cyclic->mat);
124: MatDestroy(&Zm);
125: MatDestroy(&Zn);
126: } else {
127: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&cyclic->x1);
128: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&cyclic->x2);
129: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&cyclic->y1);
130: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&cyclic->y2);
131: PetscLogObjectParent(svd,cyclic->x1);
132: PetscLogObjectParent(svd,cyclic->x2);
133: PetscLogObjectParent(svd,cyclic->y1);
134: PetscLogObjectParent(svd,cyclic->y2);
135: MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,svd,&cyclic->mat);
136: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))MatMult_Cyclic);
137: MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic);
138: }
139: PetscLogObjectParent(svd,cyclic->mat);
140: }
142: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
143: EPSSetOperators(cyclic->eps,cyclic->mat,NULL);
144: EPSSetProblemType(cyclic->eps,EPS_HEP);
145: if (svd->which == SVD_LARGEST) {
146: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
147: } else {
148: EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
149: EPSSetTarget(cyclic->eps,0.0);
150: }
151: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
152: EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
153: /* Transfer the trackall option from svd to eps */
154: SVDGetTrackAll(svd,&trackall);
155: EPSSetTrackAll(cyclic->eps,trackall);
156: /* Transfer the initial subspace from svd to eps */
157: if (svd->nini<0 || svd->ninil<0) {
158: for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
159: MatGetVecs(cyclic->mat,&v,NULL);
160: VecGetArray(v,&va);
161: if (i<-svd->ninil) {
162: VecGetSize(svd->ISL[i],&isl);
163: if (isl!=m) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
164: VecGetArrayRead(svd->ISL[i],&isa);
165: PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
166: VecRestoreArrayRead(svd->IS[i],&isa);
167: } else {
168: PetscMemzero(&va,sizeof(PetscScalar)*m);
169: }
170: if (i<-svd->nini) {
171: VecGetSize(svd->IS[i],&isl);
172: if (isl!=n) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
173: VecGetArrayRead(svd->IS[i],&isa);
174: PetscMemcpy(va+m,isa,sizeof(PetscScalar)*n);
175: VecRestoreArrayRead(svd->IS[i],&isa);
176: } else {
177: PetscMemzero(va+m,sizeof(PetscScalar)*n);
178: }
179: VecRestoreArray(v,&va);
180: VecDestroy(&svd->IS[i]);
181: svd->IS[i] = v;
182: }
183: svd->nini = PetscMin(svd->nini,svd->ninil);
184: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
185: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
186: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
187: }
188: if (cyclic->setfromoptionscalled) {
189: EPSSetFromOptions(cyclic->eps);
190: cyclic->setfromoptionscalled = PETSC_FALSE;
191: }
192: EPSSetUp(cyclic->eps);
193: EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
194: svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
195: EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);
197: if (svd->ncv != svd->n) {
198: VecDestroyVecs(svd->n,&svd->U);
199: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
200: PetscLogObjectParents(svd,svd->ncv,svd->U);
201: }
202: return(0);
203: }
207: PetscErrorCode SVDSolve_Cyclic(SVD svd)
208: {
209: PetscErrorCode ierr;
210: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
211: PetscInt i,j,M,N,m,n;
212: PetscScalar sigma;
213: const PetscScalar *px;
214: Vec x,x1,x2;
217: EPSSolve(cyclic->eps);
218: EPSGetConverged(cyclic->eps,&svd->nconv);
219: EPSGetIterationNumber(cyclic->eps,&svd->its);
220: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
222: MatGetVecs(cyclic->mat,&x,NULL);
223: SVDMatGetSize(svd,&M,&N);
224: SVDMatGetLocalSize(svd,&m,&n);
225: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&x1);
226: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&x2);
227: for (i=0,j=0;i<svd->nconv;i++) {
228: EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
229: if (PetscRealPart(sigma) > 0.0) {
230: svd->sigma[j] = PetscRealPart(sigma);
231: VecGetArrayRead(x,&px);
232: VecPlaceArray(x1,px);
233: VecPlaceArray(x2,px+m);
234: VecCopy(x1,svd->U[j]);
235: VecScale(svd->U[j],1.0/PetscSqrtReal(2.0));
236: VecCopy(x2,svd->V[j]);
237: VecScale(svd->V[j],1.0/PetscSqrtReal(2.0));
238: VecResetArray(x1);
239: VecResetArray(x2);
240: VecRestoreArrayRead(x,&px);
241: j++;
242: }
243: }
244: svd->nconv = j;
246: VecDestroy(&x);
247: VecDestroy(&x1);
248: VecDestroy(&x2);
249: return(0);
250: }
254: static PetscErrorCode SVDMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
255: {
256: PetscInt i,j;
257: SVD svd = (SVD)ctx;
258: PetscScalar er,ei;
262: nconv = 0;
263: for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
264: er = eigr[i]; ei = eigi[i];
265: STBackTransform(eps->st,1,&er,&ei);
266: if (PetscRealPart(er) > 0.0) {
267: svd->sigma[j] = PetscRealPart(er);
268: svd->errest[j] = errest[i];
269: if (errest[i] < svd->tol) nconv++;
270: j++;
271: }
272: }
273: nest = j;
274: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
275: return(0);
276: }
280: PetscErrorCode SVDSetFromOptions_Cyclic(SVD svd)
281: {
283: PetscBool set,val;
284: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
285: ST st;
288: cyclic->setfromoptionscalled = PETSC_TRUE;
289: PetscOptionsHead("SVD Cyclic Options");
290: PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
291: if (set) {
292: SVDCyclicSetExplicitMatrix(svd,val);
293: }
294: if (!cyclic->explicitmatrix) {
295: /* use as default an ST with shell matrix and Jacobi */
296: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
297: EPSGetST(cyclic->eps,&st);
298: STSetMatMode(st,ST_MATMODE_SHELL);
299: }
300: PetscOptionsTail();
301: return(0);
302: }
306: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
307: {
308: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
311: cyclic->explicitmatrix = explicitmatrix;
312: return(0);
313: }
317: /*@
318: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
319: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
321: Logically Collective on SVD
323: Input Parameters:
324: + svd - singular value solver
325: - explicit - boolean flag indicating if H(A) is built explicitly
327: Options Database Key:
328: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
330: Level: advanced
332: .seealso: SVDCyclicGetExplicitMatrix()
333: @*/
334: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
335: {
341: PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
342: return(0);
343: }
347: static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
348: {
349: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
352: *explicitmatrix = cyclic->explicitmatrix;
353: return(0);
354: }
358: /*@
359: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly
361: Not Collective
363: Input Parameter:
364: . svd - singular value solver
366: Output Parameter:
367: . explicit - the mode flag
369: Level: advanced
371: .seealso: SVDCyclicSetExplicitMatrix()
372: @*/
373: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
374: {
380: PetscTryMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
381: return(0);
382: }
386: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
387: {
388: PetscErrorCode ierr;
389: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
392: PetscObjectReference((PetscObject)eps);
393: EPSDestroy(&cyclic->eps);
394: cyclic->eps = eps;
395: PetscLogObjectParent(svd,cyclic->eps);
396: svd->setupcalled = 0;
397: return(0);
398: }
402: /*@
403: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
404: singular value solver.
406: Collective on SVD
408: Input Parameters:
409: + svd - singular value solver
410: - eps - the eigensolver object
412: Level: advanced
414: .seealso: SVDCyclicGetEPS()
415: @*/
416: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
417: {
424: PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
425: return(0);
426: }
430: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
431: {
433: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
436: if (!cyclic->eps) {
437: EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
438: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
439: EPSAppendOptionsPrefix(cyclic->eps,"svd_");
440: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
441: PetscLogObjectParent(svd,cyclic->eps);
442: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
443: EPSSetIP(cyclic->eps,svd->ip);
444: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
445: EPSMonitorSet(cyclic->eps,SVDMonitor_Cyclic,svd,NULL);
446: }
447: *eps = cyclic->eps;
448: return(0);
449: }
453: /*@
454: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
455: to the singular value solver.
457: Not Collective
459: Input Parameter:
460: . svd - singular value solver
462: Output Parameter:
463: . eps - the eigensolver object
465: Level: advanced
467: .seealso: SVDCyclicSetEPS()
468: @*/
469: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
470: {
476: PetscTryMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
477: return(0);
478: }
482: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
483: {
485: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
488: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
489: PetscViewerASCIIPrintf(viewer," Cyclic: %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
490: PetscViewerASCIIPushTab(viewer);
491: EPSView(cyclic->eps,viewer);
492: PetscViewerASCIIPopTab(viewer);
493: return(0);
494: }
498: PetscErrorCode SVDReset_Cyclic(SVD svd)
499: {
501: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
504: if (!cyclic->eps) { EPSReset(cyclic->eps); }
505: MatDestroy(&cyclic->mat);
506: VecDestroy(&cyclic->x1);
507: VecDestroy(&cyclic->x2);
508: VecDestroy(&cyclic->y1);
509: VecDestroy(&cyclic->y2);
510: return(0);
511: }
515: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
516: {
518: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
521: EPSDestroy(&cyclic->eps);
522: PetscFree(svd->data);
523: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
524: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
525: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
526: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
527: return(0);
528: }
532: PETSC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
533: {
535: SVD_CYCLIC *cyclic;
538: PetscNewLog(svd,SVD_CYCLIC,&cyclic);
539: svd->data = (void*)cyclic;
540: svd->ops->solve = SVDSolve_Cyclic;
541: svd->ops->setup = SVDSetUp_Cyclic;
542: svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
543: svd->ops->destroy = SVDDestroy_Cyclic;
544: svd->ops->reset = SVDReset_Cyclic;
545: svd->ops->view = SVDView_Cyclic;
546: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
547: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
548: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
549: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
550: return(0);
551: }