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-2007, Universidad Politecnica de Valencia, Spain
13: This file is part of SLEPc. See the README file for conditions of use
14: and additional information.
15: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16: */
18: #include src/svd/svdimpl.h
19: #include slepceps.h
21: typedef struct {
22: PetscTruth explicitmatrix;
23: EPS eps;
24: Mat mat;
25: Vec x1,x2,y1,y2;
26: } SVD_CYCLIC;
30: PetscErrorCode ShellMatMult_CYCLIC(Mat B,Vec x, Vec y)
31: {
33: SVD svd;
34: SVD_CYCLIC *cyclic;
35: PetscScalar *px,*py;
36: PetscInt m;
37:
39: MatShellGetContext(B,(void**)&svd);
40: cyclic = (SVD_CYCLIC *)svd->data;
41: SVDMatGetLocalSize(svd,&m,PETSC_NULL);
42: VecGetArray(x,&px);
43: VecGetArray(y,&py);
44: VecPlaceArray(cyclic->x1,px);
45: VecPlaceArray(cyclic->x2,px+m);
46: VecPlaceArray(cyclic->y1,py);
47: VecPlaceArray(cyclic->y2,py+m);
48: SVDMatMult(svd,PETSC_FALSE,cyclic->x2,cyclic->y1);
49: SVDMatMult(svd,PETSC_TRUE,cyclic->x1,cyclic->y2);
50: VecResetArray(cyclic->x1);
51: VecResetArray(cyclic->x2);
52: VecResetArray(cyclic->y1);
53: VecResetArray(cyclic->y2);
54: VecRestoreArray(x,&px);
55: VecRestoreArray(y,&py);
56: return(0);
57: }
61: PetscErrorCode ShellMatGetDiagonal_CYCLIC(Mat B,Vec diag)
62: {
64:
66: VecSet(diag,0.0);
67: return(0);
68: }
73: PetscErrorCode SVDSetUp_CYCLIC(SVD svd)
74: {
75: PetscErrorCode ierr;
76: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
77: PetscInt M,N,m,n,i,j,start,end,ncols,*pos;
78: const PetscInt *cols;
79: const PetscScalar *vals;
82:
83: if (cyclic->mat) {
84: MatDestroy(cyclic->mat);
85: }
86: if (cyclic->x1) {
87: VecDestroy(cyclic->x1);
88: VecDestroy(cyclic->x2);
89: VecDestroy(cyclic->y1);
90: VecDestroy(cyclic->y2);
91: }
93: SVDMatGetSize(svd,&M,&N);
94: SVDMatGetLocalSize(svd,&m,&n);
95: if (cyclic->explicitmatrix) {
96: cyclic->x1 = cyclic->x2 = cyclic->y1 = cyclic->y2 = PETSC_NULL;
97: MatCreate(svd->comm,&cyclic->mat);
98: MatSetSizes(cyclic->mat,m+n,m+n,M+N,M+N);
99: MatSetFromOptions(cyclic->mat);
100: if (svd->AT) {
101: MatGetOwnershipRange(svd->AT,&start,&end);
102: for (i=start;i<end;i++) {
103: MatGetRow(svd->AT,i,&ncols,&cols,&vals);
104: j = i + M;
105: MatSetValues(cyclic->mat,1,&j,ncols,cols,vals,INSERT_VALUES);
106: MatSetValues(cyclic->mat,ncols,cols,1,&j,vals,INSERT_VALUES);
107: MatRestoreRow(svd->AT,i,&ncols,&cols,&vals);
108: }
109: } else {
110: PetscMalloc(sizeof(PetscInt)*n,&pos);
111: MatGetOwnershipRange(svd->A,&start,&end);
112: for (i=start;i<end;i++) {
113: MatGetRow(svd->A,i,&ncols,&cols,&vals);
114: for (j=0;j<ncols;j++)
115: pos[j] = cols[j] + M;
116: MatSetValues(cyclic->mat,1,&i,ncols,pos,vals,INSERT_VALUES);
117: MatSetValues(cyclic->mat,ncols,pos,1,&i,vals,INSERT_VALUES);
118: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
119: }
120: PetscFree(pos);
121: }
122: MatAssemblyBegin(cyclic->mat,MAT_FINAL_ASSEMBLY);
123: MatAssemblyEnd(cyclic->mat,MAT_FINAL_ASSEMBLY);
124: } else {
125: VecCreateMPIWithArray(svd->comm,m,M,PETSC_NULL,&cyclic->x1);
126: VecCreateMPIWithArray(svd->comm,n,N,PETSC_NULL,&cyclic->x2);
127: VecCreateMPIWithArray(svd->comm,m,M,PETSC_NULL,&cyclic->y1);
128: VecCreateMPIWithArray(svd->comm,n,N,PETSC_NULL,&cyclic->y2);
129: MatCreateShell(svd->comm,m+n,m+n,M+N,M+N,svd,&cyclic->mat);
130: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CYCLIC);
131: MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CYCLIC);
132: }
134: EPSSetOperators(cyclic->eps,cyclic->mat,PETSC_NULL);
135: EPSSetProblemType(cyclic->eps,EPS_HEP);
136: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv);
137: EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
138: EPSSetUp(cyclic->eps);
139: EPSGetDimensions(cyclic->eps,PETSC_NULL,&svd->ncv);
140: EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);
142: if (svd->U) {
143: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
144: PetscFree(svd->U);
145: }
146: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
147: for (i=0;i<svd->ncv;i++) { SVDMatGetVecs(svd,PETSC_NULL,svd->U+i); }
149: return(0);
150: }
154: PetscErrorCode SVDSolve_CYCLIC(SVD svd)
155: {
157: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
158: int i,j;
159: PetscInt M,m,idx,start,end;
160: PetscScalar sigma,*px;
161: Vec x;
162: IS isU,isV;
163: VecScatter vsU,vsV;
164:
166: EPSSetWhichEigenpairs(cyclic->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_MAGNITUDE);
167: EPSSolve(cyclic->eps);
168: EPSGetConverged(cyclic->eps,&svd->nconv);
169: EPSGetIterationNumber(cyclic->eps,&svd->its);
170: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
172: MatGetVecs(cyclic->mat,&x,PETSC_NULL);
173: if (cyclic->explicitmatrix) {
174: EPSGetOperationCounters(cyclic->eps,&svd->matvecs,PETSC_NULL,PETSC_NULL);
175: SVDMatGetSize(svd,&M,PETSC_NULL);
176: VecGetOwnershipRange(svd->U[0],&start,&end);
177: ISCreateBlock(svd->comm,end-start,1,&start,&isU);
178: VecScatterCreate(x,isU,svd->U[0],PETSC_NULL,&vsU);
180: VecGetOwnershipRange(svd->V[0],&start,&end);
181: idx = start + M;
182: ISCreateBlock(svd->comm,end-start,1,&idx,&isV);
183: VecScatterCreate(x,isV,svd->V[0],PETSC_NULL,&vsV);
185: for (i=0,j=0;i<svd->nconv;i++) {
186: EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
187: if (PetscRealPart(sigma) > 0.0) {
188: svd->sigma[j] = PetscRealPart(sigma);
189: VecScatterBegin(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
190: VecScatterBegin(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
191: VecScatterEnd(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
192: VecScatterEnd(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
193: VecScale(svd->U[j],1.0/sqrt(2.0));
194: VecScale(svd->V[j],1.0/sqrt(2.0));
195: j++;
196: }
197: }
198:
199: ISDestroy(isU);
200: VecScatterDestroy(vsU);
201: ISDestroy(isV);
202: VecScatterDestroy(vsV);
203: } else {
204: SVDMatGetLocalSize(svd,&m,PETSC_NULL);
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: VecGetArray(x,&px);
210: VecPlaceArray(cyclic->x1,px);
211: VecPlaceArray(cyclic->x2,px+m);
212:
213: VecCopy(cyclic->x1,svd->U[j]);
214: VecScale(svd->U[j],1.0/sqrt(2.0));
216: VecCopy(cyclic->x2,svd->V[j]);
217: VecScale(svd->V[j],1.0/sqrt(2.0));
218:
219: VecResetArray(cyclic->x1);
220: VecResetArray(cyclic->x2);
221: VecRestoreArray(x,&px);
222: j++;
223: }
224: }
225: }
226: svd->nconv = j;
228: VecDestroy(x);
229: return(0);
230: }
234: PetscErrorCode SVDMonitor_CYCLIC(EPS eps,int its,int nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,int nest,void *ctx)
235: {
236: int i,j;
237: SVD svd = (SVD)ctx;
240: nconv = 0;
241: for (i=0,j=0;i<nest;i++) {
242: if (PetscRealPart(eigr[i]) > 0.0) {
243: svd->sigma[j] = PetscRealPart(eigr[i]);
244: svd->errest[j] = errest[i];
245: if (errest[i] < svd->tol) nconv++;
246: j++;
247: }
248: }
249: nest = j;
250: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
251: return(0);
252: }
256: PetscErrorCode SVDSetFromOptions_CYCLIC(SVD svd)
257: {
259: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
260: ST st;
263: PetscOptionsBegin(svd->comm,svd->prefix,"CYCLIC Singular Value Solver Options","SVD");
264: PetscOptionsTruth("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",PETSC_FALSE,&cyclic->explicitmatrix,PETSC_NULL);
265: if (cyclic->explicitmatrix) {
266: /* don't build the transpose */
267: if (svd->transmode == PETSC_DECIDE)
268: svd->transmode = SVD_TRANSPOSE_IMPLICIT;
269: } else {
270: /* use as default an ST with shell matrix and Jacobi */
271: EPSGetST(cyclic->eps,&st);
272: STSetMatMode(st,STMATMODE_SHELL);
273: }
274: PetscOptionsEnd();
275: EPSSetFromOptions(cyclic->eps);
276: PetscOptionsTail();
277: return(0);
278: }
283: PetscErrorCode SVDCyclicSetExplicitMatrix_CYCLIC(SVD svd,PetscTruth explicitmatrix)
284: {
285: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
288: cyclic->explicitmatrix = explicitmatrix;
289: return(0);
290: }
295: /*@
296: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
297: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
299: Collective on SVD
301: Input Parameters:
302: + svd - singular value solver
303: - explicit - boolean flag indicating if H(A) is built explicitly
305: Options Database Key:
306: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
308: Level: advanced
310: .seealso: SVDCyclicGetExplicitMatrix()
311: @*/
312: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscTruth explicitmatrix)
313: {
314: PetscErrorCode ierr, (*f)(SVD,PetscTruth);
318: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",(void (**)())&f);
319: if (f) {
320: (*f)(svd,explicitmatrix);
321: }
322: return(0);
323: }
328: PetscErrorCode SVDCyclicGetExplicitMatrix_CYCLIC(SVD svd,PetscTruth *explicitmatrix)
329: {
330: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
334: *explicitmatrix = cyclic->explicitmatrix;
335: return(0);
336: }
341: /*@C
342: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly
344: Not collective
346: Input Parameter:
347: . svd - singular value solver
349: Output Parameter:
350: . explicit - the mode flag
352: Level: advanced
354: .seealso: SVDCyclicSetExplicitMatrix()
355: @*/
356: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscTruth *explicitmatrix)
357: {
358: PetscErrorCode ierr, (*f)(SVD,PetscTruth*);
362: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",(void (**)())&f);
363: if (f) {
364: (*f)(svd,explicitmatrix);
365: }
366: return(0);
367: }
372: PetscErrorCode SVDCyclicSetEPS_CYCLIC(SVD svd,EPS eps)
373: {
374: PetscErrorCode ierr;
375: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
380: PetscObjectReference((PetscObject)eps);
381: EPSDestroy(cyclic->eps);
382: cyclic->eps = eps;
383: svd->setupcalled = 0;
384: return(0);
385: }
390: /*@
391: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
392: singular value solver.
394: Collective on SVD
396: Input Parameters:
397: + svd - singular value solver
398: - eps - the eigensolver object
400: Level: advanced
402: .seealso: SVDCyclicGetEPS()
403: @*/
404: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
405: {
406: PetscErrorCode ierr, (*f)(SVD,EPS eps);
410: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetEPS_C",(void (**)())&f);
411: if (f) {
412: (*f)(svd,eps);
413: }
414: return(0);
415: }
420: PetscErrorCode SVDCyclicGetEPS_CYCLIC(SVD svd,EPS *eps)
421: {
422: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
426: *eps = cyclic->eps;
427: return(0);
428: }
433: /*@C
434: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
435: to the singular value solver.
437: Not Collective
439: Input Parameter:
440: . svd - singular value solver
442: Output Parameter:
443: . eps - the eigensolver object
445: Level: advanced
447: .seealso: SVDCyclicSetEPS()
448: @*/
449: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
450: {
451: PetscErrorCode ierr, (*f)(SVD,EPS *eps);
455: PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetEPS_C",(void (**)())&f);
456: if (f) {
457: (*f)(svd,eps);
458: }
459: return(0);
460: }
464: PetscErrorCode SVDView_CYCLIC(SVD svd,PetscViewer viewer)
465: {
466: PetscErrorCode ierr;
467: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
470: if (cyclic->explicitmatrix) {
471: PetscViewerASCIIPrintf(viewer,"cyclic matrix: explicit\n");
472: } else {
473: PetscViewerASCIIPrintf(viewer,"cyclic matrix: implicit\n");
474: }
475: EPSView(cyclic->eps,viewer);
476: return(0);
477: }
481: PetscErrorCode SVDDestroy_CYCLIC(SVD svd)
482: {
483: PetscErrorCode ierr;
484: SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;
487: EPSDestroy(cyclic->eps);
488: if (cyclic->mat) { MatDestroy(cyclic->mat); }
489: if (cyclic->x1) {
490: VecDestroy(cyclic->x1);
491: VecDestroy(cyclic->x2);
492: VecDestroy(cyclic->y1);
493: VecDestroy(cyclic->y2);
494: }
495: PetscFree(svd->data);
496: return(0);
497: }
502: PetscErrorCode SVDCreate_CYCLIC(SVD svd)
503: {
504: PetscErrorCode ierr;
505: SVD_CYCLIC *cyclic;
506:
508: PetscNew(SVD_CYCLIC,&cyclic);
509: PetscLogObjectMemory(svd,sizeof(SVD_CYCLIC));
510: svd->data = (void *)cyclic;
511: svd->ops->solve = SVDSolve_CYCLIC;
512: svd->ops->setup = SVDSetUp_CYCLIC;
513: svd->ops->setfromoptions = SVDSetFromOptions_CYCLIC;
514: svd->ops->destroy = SVDDestroy_CYCLIC;
515: svd->ops->view = SVDView_CYCLIC;
516: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","SVDCyclicSetEPS_CYCLIC",SVDCyclicSetEPS_CYCLIC);
517: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","SVDCyclicGetEPS_CYCLIC",SVDCyclicGetEPS_CYCLIC);
518: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","SVDCyclicSetExplicitMatrix_CYCLIC",SVDCyclicSetExplicitMatrix_CYCLIC);
519: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","SVDCyclicGetExplicitMatrix_CYCLIC",SVDCyclicGetExplicitMatrix_CYCLIC);
521: EPSCreate(svd->comm,&cyclic->eps);
522: EPSSetOptionsPrefix(cyclic->eps,svd->prefix);
523: EPSAppendOptionsPrefix(cyclic->eps,"svd_");
524: PetscLogObjectParent(svd,cyclic->eps);
525: EPSSetIP(cyclic->eps,svd->ip);
526: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
527: EPSMonitorSet(cyclic->eps,SVDMonitor_CYCLIC,svd,PETSC_NULL);
528: cyclic->explicitmatrix = PETSC_FALSE;
529: cyclic->mat = PETSC_NULL;
530: cyclic->x1 = PETSC_NULL;
531: cyclic->x2 = PETSC_NULL;
532: cyclic->y1 = PETSC_NULL;
533: cyclic->y2 = PETSC_NULL;
534: return(0);
535: }