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-2010, Universidad 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
30: #include private/epsimpl.h
31: #include slepceps.h
33: typedef struct {
34: PetscTruth explicitmatrix;
35: EPS eps;
36: PetscTruth setfromoptionscalled;
37: Mat mat;
38: Vec x1,x2,y1,y2;
39: } SVD_CYCLIC;
43: PetscErrorCode ShellMatMult_CYCLIC(Mat B,Vec x, Vec y)
44: {
46: SVD svd;
47: SVD_CYCLIC *cyclic;
48: PetscScalar *px,*py;
49: PetscInt m;
50:
52: MatShellGetContext(B,(void**)&svd);
53: cyclic = (SVD_CYCLIC *)svd->data;
54: SVDMatGetLocalSize(svd,&m,PETSC_NULL);
55: VecGetArray(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: VecRestoreArray(x,&px);
68: VecRestoreArray(y,&py);
69: return(0);
70: }
74: PetscErrorCode ShellMatGetDiagonal_CYCLIC(Mat B,Vec diag)
75: {
77:
79: VecSet(diag,0.0);
80: return(0);
81: }
86: PetscErrorCode SVDSetUp_CYCLIC(SVD svd)
87: {
88: PetscErrorCode ierr;
89: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
90: PetscInt M,N,m,n,i,j,start,end,ncols,*pos,nloc,isl;
91: const PetscInt *cols;
92: const PetscScalar *vals;
93: PetscScalar *pU;
94: PetscTruth trackall;
95: Vec v;
96: PetscScalar *isa,*va;
99:
100: if (cyclic->mat) {
101: MatDestroy(cyclic->mat);
102: }
103: if (cyclic->x1) {
104: VecDestroy(cyclic->x1);
105: VecDestroy(cyclic->x2);
106: VecDestroy(cyclic->y1);
107: VecDestroy(cyclic->y2);
108: }
110: SVDMatGetSize(svd,&M,&N);
111: SVDMatGetLocalSize(svd,&m,&n);
112: if (cyclic->explicitmatrix) {
113: cyclic->x1 = cyclic->x2 = cyclic->y1 = cyclic->y2 = PETSC_NULL;
114: MatCreate(((PetscObject)svd)->comm,&cyclic->mat);
115: MatSetSizes(cyclic->mat,m+n,m+n,M+N,M+N);
116: MatSetFromOptions(cyclic->mat);
117: if (svd->AT) {
118: MatGetOwnershipRange(svd->AT,&start,&end);
119: for (i=start;i<end;i++) {
120: MatGetRow(svd->AT,i,&ncols,&cols,&vals);
121: j = i + M;
122: MatSetValues(cyclic->mat,1,&j,ncols,cols,vals,INSERT_VALUES);
123: MatSetValues(cyclic->mat,ncols,cols,1,&j,vals,INSERT_VALUES);
124: MatRestoreRow(svd->AT,i,&ncols,&cols,&vals);
125: }
126: } else {
127: PetscMalloc(sizeof(PetscInt)*n,&pos);
128: MatGetOwnershipRange(svd->A,&start,&end);
129: for (i=start;i<end;i++) {
130: MatGetRow(svd->A,i,&ncols,&cols,&vals);
131: for (j=0;j<ncols;j++)
132: pos[j] = cols[j] + M;
133: MatSetValues(cyclic->mat,1,&i,ncols,pos,vals,INSERT_VALUES);
134: MatSetValues(cyclic->mat,ncols,pos,1,&i,vals,INSERT_VALUES);
135: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
136: }
137: PetscFree(pos);
138: }
139: MatAssemblyBegin(cyclic->mat,MAT_FINAL_ASSEMBLY);
140: MatAssemblyEnd(cyclic->mat,MAT_FINAL_ASSEMBLY);
141: } else {
142: VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->x1);
143: VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->x2);
144: VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->y1);
145: VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->y2);
146: MatCreateShell(((PetscObject)svd)->comm,m+n,m+n,M+N,M+N,svd,&cyclic->mat);
147: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CYCLIC);
148: MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CYCLIC);
149: }
151: EPSSetOperators(cyclic->eps,cyclic->mat,PETSC_NULL);
152: EPSSetProblemType(cyclic->eps,EPS_HEP);
153: EPSSetWhichEigenpairs(cyclic->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_MAGNITUDE);
154: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
155: EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
156: /* Transfer the trackall option from svd to eps */
157: SVDGetTrackAll(svd,&trackall);
158: EPSSetTrackAll(cyclic->eps,trackall);
159: /* Transfer the initial subspace from svd to eps */
160: if (svd->nini < 0) {
161: for (i=0; i<-svd->nini; i++) {
162: MatGetVecs(cyclic->mat,&v,PETSC_NULL);
163: VecGetArray(v,&va);
164: VecGetArray(svd->IS[i],&isa);
165: VecGetSize(svd->IS[i],&isl);
166: if (isl == m) {
167: PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
168: PetscMemzero(&va[m],sizeof(PetscScalar)*n);
169: } else if (isl == n) {
170: PetscMemzero(va,sizeof(PetscScalar)*m);
171: PetscMemcpy(&va[m],isa,sizeof(PetscScalar)*n);
172: } else {
173: SETERRQ(PETSC_ERR_SUP,"Size of the initial subspace vectors should match to some dimension of A");
174: }
175: VecRestoreArray(v,&va);
176: VecRestoreArray(svd->IS[i],&isa);
177: VecDestroy(svd->IS[i]);
178: svd->IS[i] = v;
179: }
180: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
181: for (i=0; i<-svd->nini; i++) {
182: VecDestroy(svd->IS[i]);
183: }
184: PetscFree(svd->IS);
185: svd->IS = PETSC_NULL;
186: svd->nini = 0;
187: }
188: if (cyclic->setfromoptionscalled) {
189: EPSSetFromOptions(cyclic->eps);
190: cyclic->setfromoptionscalled = PETSC_FALSE;
191: }
192: EPSSetUp(cyclic->eps);
193: EPSGetDimensions(cyclic->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
194: EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);
196: if (svd->ncv != svd->n) {
197: if (svd->U) {
198: VecGetArray(svd->U[0],&pU);
199: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
200: PetscFree(pU);
201: PetscFree(svd->U);
202: }
203: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
204: SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
205: PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
206: for (i=0;i<svd->ncv;i++) {
207: VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
208: }
209: }
211: return(0);
212: }
216: PetscErrorCode SVDSolve_CYCLIC(SVD svd)
217: {
219: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
220: PetscInt i,j,M,m,idx,start,end;
221: PetscScalar sigma,*px;
222: Vec x;
223: IS isU,isV;
224: VecScatter vsU,vsV;
225:
227: EPSSolve(cyclic->eps);
228: EPSGetConverged(cyclic->eps,&svd->nconv);
229: EPSGetIterationNumber(cyclic->eps,&svd->its);
230: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
232: MatGetVecs(cyclic->mat,&x,PETSC_NULL);
233: if (cyclic->explicitmatrix) {
234: EPSGetOperationCounters(cyclic->eps,&svd->matvecs,PETSC_NULL,PETSC_NULL);
235: SVDMatGetSize(svd,&M,PETSC_NULL);
236: VecGetOwnershipRange(svd->U[0],&start,&end);
237: ISCreateBlock(((PetscObject)svd)->comm,end-start,1,&start,&isU);
238: VecScatterCreate(x,isU,svd->U[0],PETSC_NULL,&vsU);
240: VecGetOwnershipRange(svd->V[0],&start,&end);
241: idx = start + M;
242: ISCreateBlock(((PetscObject)svd)->comm,end-start,1,&idx,&isV);
243: VecScatterCreate(x,isV,svd->V[0],PETSC_NULL,&vsV);
245: for (i=0,j=0;i<svd->nconv;i++) {
246: EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
247: if (PetscRealPart(sigma) > 0.0) {
248: svd->sigma[j] = PetscRealPart(sigma);
249: VecScatterBegin(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
250: VecScatterBegin(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
251: VecScatterEnd(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
252: VecScatterEnd(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
253: VecScale(svd->U[j],1.0/sqrt(2.0));
254: VecScale(svd->V[j],1.0/sqrt(2.0));
255: j++;
256: }
257: }
258:
259: ISDestroy(isU);
260: VecScatterDestroy(vsU);
261: ISDestroy(isV);
262: VecScatterDestroy(vsV);
263: } else {
264: SVDMatGetLocalSize(svd,&m,PETSC_NULL);
265: for (i=0,j=0;i<svd->nconv;i++) {
266: EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
267: if (PetscRealPart(sigma) > 0.0) {
268: svd->sigma[j] = PetscRealPart(sigma);
269: VecGetArray(x,&px);
270: VecPlaceArray(cyclic->x1,px);
271: VecPlaceArray(cyclic->x2,px+m);
272:
273: VecCopy(cyclic->x1,svd->U[j]);
274: VecScale(svd->U[j],1.0/sqrt(2.0));
276: VecCopy(cyclic->x2,svd->V[j]);
277: VecScale(svd->V[j],1.0/sqrt(2.0));
278:
279: VecResetArray(cyclic->x1);
280: VecResetArray(cyclic->x2);
281: VecRestoreArray(x,&px);
282: j++;
283: }
284: }
285: }
286: svd->nconv = j;
288: VecDestroy(x);
289: return(0);
290: }
294: PetscErrorCode SVDMonitor_CYCLIC(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
295: {
296: PetscInt i,j;
297: SVD svd = (SVD)ctx;
298: PetscScalar er,ei;
302: nconv = 0;
303: for (i=0,j=0;i<nest;i++) {
304: er = eigr[i]; ei = eigi[i];
305: STBackTransform(eps->OP, 1, &er, &ei);
306: if (PetscRealPart(er) > 0.0) {
307: svd->sigma[j] = PetscRealPart(er);
308: svd->errest[j] = errest[i];
309: if (errest[i] < svd->tol) nconv++;
310: j++;
311: }
312: }
313: nest = j;
314: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
315: return(0);
316: }
320: PetscErrorCode SVDSetFromOptions_CYCLIC(SVD svd)
321: {
323: PetscTruth set,val;
324: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
325: ST st;
328: PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"CYCLIC Singular Value Solver Options","SVD");
329: PetscOptionsTruth("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
330: if (set) {
331: SVDCyclicSetExplicitMatrix(svd,val);
332: }
333: if (cyclic->explicitmatrix) {
334: /* don't build the transpose */
335: if (svd->transmode == PETSC_DECIDE)
336: svd->transmode = SVD_TRANSPOSE_IMPLICIT;
337: } else {
338: /* use as default an ST with shell matrix and Jacobi */
339: EPSGetST(cyclic->eps,&st);
340: STSetMatMode(st,ST_MATMODE_SHELL);
341: }
342: PetscOptionsEnd();
343: cyclic->setfromoptionscalled = PETSC_TRUE;
344: return(0);
345: }
350: PetscErrorCode SVDCyclicSetExplicitMatrix_CYCLIC(SVD svd,PetscTruth explicitmatrix)
351: {
352: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
355: cyclic->explicitmatrix = explicitmatrix;
356: return(0);
357: }
362: /*@
363: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
364: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
366: Collective on SVD
368: Input Parameters:
369: + svd - singular value solver
370: - explicit - boolean flag indicating if H(A) is built explicitly
372: Options Database Key:
373: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
375: Level: advanced
377: .seealso: SVDCyclicGetExplicitMatrix()
378: @*/
379: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscTruth explicitmatrix)
380: {
381: PetscErrorCode ierr, (*f)(SVD,PetscTruth);
385: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",(void (**)())&f);
386: if (f) {
387: (*f)(svd,explicitmatrix);
388: }
389: return(0);
390: }
395: PetscErrorCode SVDCyclicGetExplicitMatrix_CYCLIC(SVD svd,PetscTruth *explicitmatrix)
396: {
397: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
401: *explicitmatrix = cyclic->explicitmatrix;
402: return(0);
403: }
408: /*@
409: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly
411: Not collective
413: Input Parameter:
414: . svd - singular value solver
416: Output Parameter:
417: . explicit - the mode flag
419: Level: advanced
421: .seealso: SVDCyclicSetExplicitMatrix()
422: @*/
423: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscTruth *explicitmatrix)
424: {
425: PetscErrorCode ierr, (*f)(SVD,PetscTruth*);
429: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",(void (**)())&f);
430: if (f) {
431: (*f)(svd,explicitmatrix);
432: }
433: return(0);
434: }
439: PetscErrorCode SVDCyclicSetEPS_CYCLIC(SVD svd,EPS eps)
440: {
441: PetscErrorCode ierr;
442: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
447: PetscObjectReference((PetscObject)eps);
448: EPSDestroy(cyclic->eps);
449: cyclic->eps = eps;
450: svd->setupcalled = 0;
451: return(0);
452: }
457: /*@
458: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
459: singular value solver.
461: Collective on SVD
463: Input Parameters:
464: + svd - singular value solver
465: - eps - the eigensolver object
467: Level: advanced
469: .seealso: SVDCyclicGetEPS()
470: @*/
471: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
472: {
473: PetscErrorCode ierr, (*f)(SVD,EPS eps);
477: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetEPS_C",(void (**)())&f);
478: if (f) {
479: (*f)(svd,eps);
480: }
481: return(0);
482: }
487: PetscErrorCode SVDCyclicGetEPS_CYCLIC(SVD svd,EPS *eps)
488: {
489: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
493: *eps = cyclic->eps;
494: return(0);
495: }
500: /*@
501: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
502: to the singular value solver.
504: Not Collective
506: Input Parameter:
507: . svd - singular value solver
509: Output Parameter:
510: . eps - the eigensolver object
512: Level: advanced
514: .seealso: SVDCyclicSetEPS()
515: @*/
516: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
517: {
518: PetscErrorCode ierr, (*f)(SVD,EPS *eps);
522: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetEPS_C",(void (**)())&f);
523: if (f) {
524: (*f)(svd,eps);
525: }
526: return(0);
527: }
531: PetscErrorCode SVDView_CYCLIC(SVD svd,PetscViewer viewer)
532: {
533: PetscErrorCode ierr;
534: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
537: if (cyclic->explicitmatrix) {
538: PetscViewerASCIIPrintf(viewer,"cyclic matrix: explicit\n");
539: } else {
540: PetscViewerASCIIPrintf(viewer,"cyclic matrix: implicit\n");
541: }
542: EPSView(cyclic->eps,viewer);
543: return(0);
544: }
548: PetscErrorCode SVDDestroy_CYCLIC(SVD svd)
549: {
550: PetscErrorCode ierr;
551: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
554: EPSDestroy(cyclic->eps);
555: if (cyclic->mat) { MatDestroy(cyclic->mat); }
556: if (cyclic->x1) {
557: VecDestroy(cyclic->x1);
558: VecDestroy(cyclic->x2);
559: VecDestroy(cyclic->y1);
560: VecDestroy(cyclic->y2);
561: }
562: PetscFree(svd->data);
563: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","",PETSC_NULL);
564: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","",PETSC_NULL);
565: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","",PETSC_NULL);
566: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","",PETSC_NULL);
567: return(0);
568: }
573: PetscErrorCode SVDCreate_CYCLIC(SVD svd)
574: {
575: PetscErrorCode ierr;
576: SVD_CYCLIC *cyclic;
577:
579: PetscNew(SVD_CYCLIC,&cyclic);
580: PetscLogObjectMemory(svd,sizeof(SVD_CYCLIC));
581: svd->data = (void *)cyclic;
582: svd->ops->solve = SVDSolve_CYCLIC;
583: svd->ops->setup = SVDSetUp_CYCLIC;
584: svd->ops->setfromoptions = SVDSetFromOptions_CYCLIC;
585: svd->ops->destroy = SVDDestroy_CYCLIC;
586: svd->ops->view = SVDView_CYCLIC;
587: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","SVDCyclicSetEPS_CYCLIC",SVDCyclicSetEPS_CYCLIC);
588: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","SVDCyclicGetEPS_CYCLIC",SVDCyclicGetEPS_CYCLIC);
589: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","SVDCyclicSetExplicitMatrix_CYCLIC",SVDCyclicSetExplicitMatrix_CYCLIC);
590: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","SVDCyclicGetExplicitMatrix_CYCLIC",SVDCyclicGetExplicitMatrix_CYCLIC);
592: EPSCreate(((PetscObject)svd)->comm,&cyclic->eps);
593: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
594: EPSAppendOptionsPrefix(cyclic->eps,"svd_");
595: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
596: PetscLogObjectParent(svd,cyclic->eps);
597: EPSSetIP(cyclic->eps,svd->ip);
598: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
599: EPSMonitorSet(cyclic->eps,SVDMonitor_CYCLIC,svd,PETSC_NULL);
600: cyclic->explicitmatrix = PETSC_FALSE;
601: cyclic->mat = PETSC_NULL;
602: cyclic->x1 = PETSC_NULL;
603: cyclic->x2 = PETSC_NULL;
604: cyclic->y1 = PETSC_NULL;
605: cyclic->y2 = PETSC_NULL;
606: cyclic->setfromoptionscalled = PETSC_FALSE;
607: return(0);
608: }