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-2012, 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 <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 ShellMatMult_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;
50:
52: MatShellGetContext(B,(void**)&svd);
53: cyclic = (SVD_CYCLIC*)svd->data;
54: SVDMatGetLocalSize(svd,&m,PETSC_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 ShellMatGetDiagonal_Cyclic(Mat B,Vec diag)
75: {
77:
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;
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(((PetscObject)svd)->comm,PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
102: MatCreate(((PetscObject)svd)->comm,&Zm);
103: MatSetSizes(Zm,m,m,M,M);
104: MatSetFromOptions(Zm);
105: MatSetUp(Zm);
106: MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
107: MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
108: MatCreate(((PetscObject)svd)->comm,&Zn);
109: MatSetSizes(Zn,n,n,N,N);
110: MatSetFromOptions(Zn);
111: MatSetUp(Zn);
112: MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
114: SlepcMatTile(0.0,Zm,1.0,svd->A,1.0,svd->AT,0.0,Zn,&cyclic->mat);
115: MatDestroy(&Zm);
116: MatDestroy(&Zn);
117: } else {
118: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,m,M,PETSC_NULL,&cyclic->x1);
119: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,n,N,PETSC_NULL,&cyclic->x2);
120: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,m,M,PETSC_NULL,&cyclic->y1);
121: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,n,N,PETSC_NULL,&cyclic->y2);
122: MatCreateShell(((PetscObject)svd)->comm,m+n,m+n,M+N,M+N,svd,&cyclic->mat);
123: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))ShellMatMult_Cyclic);
124: MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_Cyclic);
125: }
126: }
128: EPSSetOperators(cyclic->eps,cyclic->mat,PETSC_NULL);
129: EPSSetProblemType(cyclic->eps,EPS_HEP);
130: if (svd->which == SVD_LARGEST) {
131: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
132: } else {
133: EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPositiveReal,PETSC_NULL);
134: EPSSetTarget(cyclic->eps,0.0);
135: }
136: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
137: EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
138: /* Transfer the trackall option from svd to eps */
139: SVDGetTrackAll(svd,&trackall);
140: EPSSetTrackAll(cyclic->eps,trackall);
141: /* Transfer the initial subspace from svd to eps */
142: if (svd->nini < 0) {
143: for (i=0; i<-svd->nini; i++) {
144: MatGetVecs(cyclic->mat,&v,PETSC_NULL);
145: VecGetArray(v,&va);
146: VecGetArrayRead(svd->IS[i],&isa);
147: VecGetSize(svd->IS[i],&isl);
148: if (isl == m) {
149: PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
150: PetscMemzero(&va[m],sizeof(PetscScalar)*n);
151: } else if (isl == n) {
152: PetscMemzero(va,sizeof(PetscScalar)*m);
153: PetscMemcpy(&va[m],isa,sizeof(PetscScalar)*n);
154: } else SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_SUP,"Size of the initial subspace vectors should match to some dimension of A");
155: VecRestoreArray(v,&va);
156: VecRestoreArrayRead(svd->IS[i],&isa);
157: VecDestroy(&svd->IS[i]);
158: svd->IS[i] = v;
159: }
160: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
161: for (i=0; i<-svd->nini; i++) {
162: VecDestroy(&svd->IS[i]);
163: }
164: PetscFree(svd->IS);
165: svd->nini = 0;
166: }
167: if (cyclic->setfromoptionscalled) {
168: EPSSetFromOptions(cyclic->eps);
169: cyclic->setfromoptionscalled = PETSC_FALSE;
170: }
171: EPSSetUp(cyclic->eps);
172: EPSGetDimensions(cyclic->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
173: svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
174: EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);
176: if (svd->ncv != svd->n) {
177: VecDestroyVecs(svd->n,&svd->U);
178: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
179: }
180: return(0);
181: }
185: PetscErrorCode SVDSolve_Cyclic(SVD svd)
186: {
187: PetscErrorCode ierr;
188: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
189: PetscInt i,j,M,N,m,n;
190: PetscScalar sigma;
191: const PetscScalar *px;
192: Vec x,x1,x2;
193:
195: EPSSolve(cyclic->eps);
196: EPSGetConverged(cyclic->eps,&svd->nconv);
197: EPSGetIterationNumber(cyclic->eps,&svd->its);
198: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
200: MatGetVecs(cyclic->mat,&x,PETSC_NULL);
201: SVDMatGetSize(svd,&M,&N);
202: SVDMatGetLocalSize(svd,&m,&n);
203: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,m,M,PETSC_NULL,&x1);
204: VecCreateMPIWithArray(((PetscObject)svd)->comm,1,n,N,PETSC_NULL,&x2);
205: for (i=0,j=0;i<svd->nconv;i++) {
206: EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
207: if (PetscRealPart(sigma) > 0.0) {
208: svd->sigma[j] = PetscRealPart(sigma);
209: VecGetArrayRead(x,&px);
210: VecPlaceArray(x1,px);
211: VecPlaceArray(x2,px+m);
212: VecCopy(x1,svd->U[j]);
213: VecScale(svd->U[j],1.0/PetscSqrtReal(2.0));
214: VecCopy(x2,svd->V[j]);
215: VecScale(svd->V[j],1.0/PetscSqrtReal(2.0));
216: VecResetArray(x1);
217: VecResetArray(x2);
218: VecRestoreArrayRead(x,&px);
219: j++;
220: }
221: }
222: svd->nconv = j;
224: VecDestroy(&x);
225: VecDestroy(&x1);
226: VecDestroy(&x2);
227: return(0);
228: }
232: static PetscErrorCode SVDMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
233: {
234: PetscInt i,j;
235: SVD svd = (SVD)ctx;
236: PetscScalar er,ei;
240: nconv = 0;
241: for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
242: er = eigr[i]; ei = eigi[i];
243: STBackTransform(eps->OP,1,&er,&ei);
244: if (PetscRealPart(er) > 0.0) {
245: svd->sigma[j] = PetscRealPart(er);
246: svd->errest[j] = errest[i];
247: if (errest[i] < svd->tol) nconv++;
248: j++;
249: }
250: }
251: nest = j;
252: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
253: return(0);
254: }
258: PetscErrorCode SVDSetFromOptions_Cyclic(SVD svd)
259: {
261: PetscBool set,val;
262: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
263: ST st;
266: cyclic->setfromoptionscalled = PETSC_TRUE;
267: PetscOptionsHead("SVD Cyclic Options");
268: PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
269: if (set) {
270: SVDCyclicSetExplicitMatrix(svd,val);
271: }
272: if (!cyclic->explicitmatrix) {
273: /* use as default an ST with shell matrix and Jacobi */
274: EPSGetST(cyclic->eps,&st);
275: STSetMatMode(st,ST_MATMODE_SHELL);
276: }
277: PetscOptionsTail();
278: return(0);
279: }
281: EXTERN_C_BEGIN
284: PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
285: {
286: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
289: cyclic->explicitmatrix = explicitmatrix;
290: return(0);
291: }
292: EXTERN_C_END
296: /*@
297: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
298: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
300: Logically Collective on SVD
302: Input Parameters:
303: + svd - singular value solver
304: - explicit - boolean flag indicating if H(A) is built explicitly
306: Options Database Key:
307: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
309: Level: advanced
311: .seealso: SVDCyclicGetExplicitMatrix()
312: @*/
313: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
314: {
320: PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
321: return(0);
322: }
324: EXTERN_C_BEGIN
327: PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
328: {
329: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
332: *explicitmatrix = cyclic->explicitmatrix;
333: return(0);
334: }
335: EXTERN_C_END
339: /*@
340: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly
342: Not Collective
344: Input Parameter:
345: . svd - singular value solver
347: Output Parameter:
348: . explicit - the mode flag
350: Level: advanced
352: .seealso: SVDCyclicSetExplicitMatrix()
353: @*/
354: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
355: {
361: PetscTryMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
362: return(0);
363: }
365: EXTERN_C_BEGIN
368: PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
369: {
370: PetscErrorCode ierr;
371: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
374: PetscObjectReference((PetscObject)eps);
375: EPSDestroy(&cyclic->eps);
376: cyclic->eps = eps;
377: PetscLogObjectParent(svd,cyclic->eps);
378: svd->setupcalled = 0;
379: return(0);
380: }
381: EXTERN_C_END
385: /*@
386: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
387: singular value solver.
389: Collective on SVD
391: Input Parameters:
392: + svd - singular value solver
393: - eps - the eigensolver object
395: Level: advanced
397: .seealso: SVDCyclicGetEPS()
398: @*/
399: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
400: {
407: PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
408: return(0);
409: }
411: EXTERN_C_BEGIN
414: PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
415: {
416: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
419: *eps = cyclic->eps;
420: return(0);
421: }
422: EXTERN_C_END
426: /*@
427: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
428: to the singular value solver.
430: Not Collective
432: Input Parameter:
433: . svd - singular value solver
435: Output Parameter:
436: . eps - the eigensolver object
438: Level: advanced
440: .seealso: SVDCyclicSetEPS()
441: @*/
442: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
443: {
449: PetscTryMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
450: return(0);
451: }
455: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
456: {
458: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
461: PetscViewerASCIIPrintf(viewer," Cyclic: %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
462: PetscViewerASCIIPushTab(viewer);
463: EPSView(cyclic->eps,viewer);
464: PetscViewerASCIIPopTab(viewer);
465: return(0);
466: }
470: PetscErrorCode SVDReset_Cyclic(SVD svd)
471: {
473: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
476: EPSReset(cyclic->eps);
477: MatDestroy(&cyclic->mat);
478: VecDestroy(&cyclic->x1);
479: VecDestroy(&cyclic->x2);
480: VecDestroy(&cyclic->y1);
481: VecDestroy(&cyclic->y2);
482: return(0);
483: }
487: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
488: {
490: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
493: EPSDestroy(&cyclic->eps);
494: PetscFree(svd->data);
495: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","",PETSC_NULL);
496: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","",PETSC_NULL);
497: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","",PETSC_NULL);
498: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","",PETSC_NULL);
499: return(0);
500: }
502: EXTERN_C_BEGIN
505: PetscErrorCode SVDCreate_Cyclic(SVD svd)
506: {
508: SVD_CYCLIC *cyclic;
509:
511: PetscNewLog(svd,SVD_CYCLIC,&cyclic);
512: svd->data = (void *)cyclic;
513: svd->ops->solve = SVDSolve_Cyclic;
514: svd->ops->setup = SVDSetUp_Cyclic;
515: svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
516: svd->ops->destroy = SVDDestroy_Cyclic;
517: svd->ops->reset = SVDReset_Cyclic;
518: svd->ops->view = SVDView_Cyclic;
519: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","SVDCyclicSetEPS_Cyclic",SVDCyclicSetEPS_Cyclic);
520: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","SVDCyclicGetEPS_Cyclic",SVDCyclicGetEPS_Cyclic);
521: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","SVDCyclicSetExplicitMatrix_Cyclic",SVDCyclicSetExplicitMatrix_Cyclic);
522: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","SVDCyclicGetExplicitMatrix_Cyclic",SVDCyclicGetExplicitMatrix_Cyclic);
524: EPSCreate(((PetscObject)svd)->comm,&cyclic->eps);
525: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
526: EPSAppendOptionsPrefix(cyclic->eps,"svd_");
527: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
528: PetscLogObjectParent(svd,cyclic->eps);
529: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
530: EPSSetIP(cyclic->eps,svd->ip);
531: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
532: EPSMonitorSet(cyclic->eps,SVDMonitor_Cyclic,svd,PETSC_NULL);
533: cyclic->explicitmatrix = PETSC_FALSE;
534: cyclic->mat = PETSC_NULL;
535: cyclic->x1 = PETSC_NULL;
536: cyclic->x2 = PETSC_NULL;
537: cyclic->y1 = PETSC_NULL;
538: cyclic->y2 = PETSC_NULL;
539: cyclic->setfromoptionscalled = PETSC_FALSE;
540: return(0);
541: }
542: EXTERN_C_END