Actual source code: cyclic.c
slepc-3.7.2 2016-07-19
1: /*
3: SLEPc singular value solver: "cyclic"
5: Method: Uses a Hermitian eigensolver for H(A) = [ 0 A ; A^T 0 ]
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2016, 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: PetscBool explicitmatrix;
32: EPS eps;
33: Mat mat;
34: Vec x1,x2,y1,y2;
35: } SVD_CYCLIC;
39: static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
40: {
41: PetscErrorCode ierr;
42: SVD svd;
43: SVD_CYCLIC *cyclic;
44: const PetscScalar *px;
45: PetscScalar *py;
46: PetscInt m;
49: MatShellGetContext(B,(void**)&svd);
50: cyclic = (SVD_CYCLIC*)svd->data;
51: SVDMatGetLocalSize(svd,&m,NULL);
52: VecGetArrayRead(x,&px);
53: VecGetArray(y,&py);
54: VecPlaceArray(cyclic->x1,px);
55: VecPlaceArray(cyclic->x2,px+m);
56: VecPlaceArray(cyclic->y1,py);
57: VecPlaceArray(cyclic->y2,py+m);
58: SVDMatMult(svd,PETSC_FALSE,cyclic->x2,cyclic->y1);
59: SVDMatMult(svd,PETSC_TRUE,cyclic->x1,cyclic->y2);
60: VecResetArray(cyclic->x1);
61: VecResetArray(cyclic->x2);
62: VecResetArray(cyclic->y1);
63: VecResetArray(cyclic->y2);
64: VecRestoreArrayRead(x,&px);
65: VecRestoreArray(y,&py);
66: return(0);
67: }
71: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
72: {
76: VecSet(diag,0.0);
77: return(0);
78: }
82: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
83: {
84: PetscErrorCode ierr;
85: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
86: PetscInt M,N,m,n,i,isl,Istart,Iend;
87: const PetscScalar *isa;
88: PetscScalar *va;
89: PetscBool trackall,gpu,issinv;
90: Vec v;
91: Mat Zm,Zn;
92: ST st;
95: PetscObjectTypeCompareAny((PetscObject)svd->A,&gpu,MATSEQAIJCUSP,MATMPIAIJCUSP,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
96: if (gpu) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Solver not implemented for GPU matrices");
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((PetscObject)svd,(PetscObject)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((PetscObject)svd,(PetscObject)cyclic->x1);
132: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->x2);
133: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->y1);
134: PetscLogObjectParent((PetscObject)svd,(PetscObject)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((PetscObject)svd,(PetscObject)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: EPSGetST(cyclic->eps,&st);
147: PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv);
148: if (issinv) {
149: EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE);
150: } else {
151: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
152: }
153: } else {
154: EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
155: EPSSetTarget(cyclic->eps,0.0);
156: }
157: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv?svd->ncv:PETSC_DEFAULT,svd->mpd?svd->mpd:PETSC_DEFAULT);
158: EPSSetTolerances(cyclic->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it?svd->max_it:PETSC_DEFAULT);
159: switch (svd->conv) {
160: case SVD_CONV_ABS:
161: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_ABS);break;
162: case SVD_CONV_REL:
163: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_REL);break;
164: case SVD_CONV_USER:
165: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
166: }
167: if (svd->stop!=SVD_STOP_BASIC) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined stopping test not supported in this solver");
168: /* Transfer the trackall option from svd to eps */
169: SVDGetTrackAll(svd,&trackall);
170: EPSSetTrackAll(cyclic->eps,trackall);
171: /* Transfer the initial subspace from svd to eps */
172: if (svd->nini<0 || svd->ninil<0) {
173: for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
174: MatCreateVecs(cyclic->mat,&v,NULL);
175: VecGetArray(v,&va);
176: if (i<-svd->ninil) {
177: VecGetSize(svd->ISL[i],&isl);
178: if (isl!=m) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
179: VecGetArrayRead(svd->ISL[i],&isa);
180: PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
181: VecRestoreArrayRead(svd->IS[i],&isa);
182: } else {
183: PetscMemzero(&va,sizeof(PetscScalar)*m);
184: }
185: if (i<-svd->nini) {
186: VecGetSize(svd->IS[i],&isl);
187: if (isl!=n) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
188: VecGetArrayRead(svd->IS[i],&isa);
189: PetscMemcpy(va+m,isa,sizeof(PetscScalar)*n);
190: VecRestoreArrayRead(svd->IS[i],&isa);
191: } else {
192: PetscMemzero(va+m,sizeof(PetscScalar)*n);
193: }
194: VecRestoreArray(v,&va);
195: VecDestroy(&svd->IS[i]);
196: svd->IS[i] = v;
197: }
198: svd->nini = PetscMin(svd->nini,svd->ninil);
199: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
200: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
201: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
202: }
203: EPSSetUp(cyclic->eps);
204: EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
205: svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
206: EPSGetTolerances(cyclic->eps,NULL,&svd->max_it);
207: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
209: svd->leftbasis = PETSC_TRUE;
210: SVDAllocateSolution(svd,0);
211: return(0);
212: }
216: PetscErrorCode SVDSolve_Cyclic(SVD svd)
217: {
218: PetscErrorCode ierr;
219: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
220: PetscInt i,j,M,N,m,n;
221: PetscScalar sigma;
222: const PetscScalar *px;
223: Vec x,x1,x2;
226: EPSSolve(cyclic->eps);
227: EPSGetConverged(cyclic->eps,&svd->nconv);
228: EPSGetIterationNumber(cyclic->eps,&svd->its);
229: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
231: MatCreateVecs(cyclic->mat,&x,NULL);
232: SVDMatGetSize(svd,&M,&N);
233: SVDMatGetLocalSize(svd,&m,&n);
234: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&x1);
235: VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&x2);
236: for (i=0,j=0;i<svd->nconv;i++) {
237: EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
238: if (PetscRealPart(sigma) > 0.0) {
239: svd->sigma[j] = PetscRealPart(sigma);
240: VecGetArrayRead(x,&px);
241: VecPlaceArray(x1,px);
242: VecPlaceArray(x2,px+m);
243: BVInsertVec(svd->U,j,x1);
244: BVScaleColumn(svd->U,j,1.0/PetscSqrtReal(2.0));
245: BVInsertVec(svd->V,j,x2);
246: BVScaleColumn(svd->V,j,1.0/PetscSqrtReal(2.0));
247: VecResetArray(x1);
248: VecResetArray(x2);
249: VecRestoreArrayRead(x,&px);
250: j++;
251: }
252: }
253: svd->nconv = j;
255: VecDestroy(&x);
256: VecDestroy(&x1);
257: VecDestroy(&x2);
258: return(0);
259: }
263: static PetscErrorCode EPSMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
264: {
265: PetscInt i,j;
266: SVD svd = (SVD)ctx;
267: PetscScalar er,ei;
271: nconv = 0;
272: for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
273: er = eigr[i]; ei = eigi[i];
274: STBackTransform(eps->st,1,&er,&ei);
275: if (PetscRealPart(er) > 0.0) {
276: svd->sigma[j] = PetscRealPart(er);
277: svd->errest[j] = errest[i];
278: if (errest[i] && errest[i] < svd->tol) nconv++;
279: j++;
280: }
281: }
282: nest = j;
283: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
284: return(0);
285: }
289: PetscErrorCode SVDSetFromOptions_Cyclic(PetscOptionItems *PetscOptionsObject,SVD svd)
290: {
292: PetscBool set,val;
293: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
294: ST st;
297: PetscOptionsHead(PetscOptionsObject,"SVD Cyclic Options");
298: PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
299: if (set) {
300: SVDCyclicSetExplicitMatrix(svd,val);
301: }
302: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
303: EPSSetFromOptions(cyclic->eps);
304: if (!cyclic->explicitmatrix) {
305: /* use as default an ST with shell matrix and Jacobi */
306: EPSGetST(cyclic->eps,&st);
307: STSetMatMode(st,ST_MATMODE_SHELL);
308: }
309: PetscOptionsTail();
310: return(0);
311: }
315: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
316: {
317: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
320: cyclic->explicitmatrix = explicitmatrix;
321: return(0);
322: }
326: /*@
327: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
328: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
330: Logically Collective on SVD
332: Input Parameters:
333: + svd - singular value solver
334: - explicit - boolean flag indicating if H(A) is built explicitly
336: Options Database Key:
337: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
339: Level: advanced
341: .seealso: SVDCyclicGetExplicitMatrix()
342: @*/
343: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
344: {
350: PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
351: return(0);
352: }
356: static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
357: {
358: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
361: *explicitmatrix = cyclic->explicitmatrix;
362: return(0);
363: }
367: /*@
368: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly
370: Not Collective
372: Input Parameter:
373: . svd - singular value solver
375: Output Parameter:
376: . explicit - the mode flag
378: Level: advanced
380: .seealso: SVDCyclicSetExplicitMatrix()
381: @*/
382: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
383: {
389: PetscUseMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
390: return(0);
391: }
395: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
396: {
397: PetscErrorCode ierr;
398: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
401: PetscObjectReference((PetscObject)eps);
402: EPSDestroy(&cyclic->eps);
403: cyclic->eps = eps;
404: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
405: svd->state = SVD_STATE_INITIAL;
406: return(0);
407: }
411: /*@
412: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
413: singular value solver.
415: Collective on SVD
417: Input Parameters:
418: + svd - singular value solver
419: - eps - the eigensolver object
421: Level: advanced
423: .seealso: SVDCyclicGetEPS()
424: @*/
425: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
426: {
433: PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
434: return(0);
435: }
439: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
440: {
442: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
445: if (!cyclic->eps) {
446: EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
447: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
448: EPSAppendOptionsPrefix(cyclic->eps,"svd_cyclic_");
449: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
450: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
451: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
452: EPSMonitorSet(cyclic->eps,EPSMonitor_Cyclic,svd,NULL);
453: }
454: *eps = cyclic->eps;
455: return(0);
456: }
460: /*@
461: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
462: to the singular value solver.
464: Not Collective
466: Input Parameter:
467: . svd - singular value solver
469: Output Parameter:
470: . eps - the eigensolver object
472: Level: advanced
474: .seealso: SVDCyclicSetEPS()
475: @*/
476: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
477: {
483: PetscUseMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
484: return(0);
485: }
489: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
490: {
492: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
493: PetscBool isascii;
496: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
497: if (isascii) {
498: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
499: PetscViewerASCIIPrintf(viewer," Cyclic: %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
500: PetscViewerASCIIPushTab(viewer);
501: EPSView(cyclic->eps,viewer);
502: PetscViewerASCIIPopTab(viewer);
503: }
504: return(0);
505: }
509: PetscErrorCode SVDReset_Cyclic(SVD svd)
510: {
512: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
515: if (!cyclic->eps) { EPSReset(cyclic->eps); }
516: MatDestroy(&cyclic->mat);
517: VecDestroy(&cyclic->x1);
518: VecDestroy(&cyclic->x2);
519: VecDestroy(&cyclic->y1);
520: VecDestroy(&cyclic->y2);
521: return(0);
522: }
526: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
527: {
529: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
532: EPSDestroy(&cyclic->eps);
533: PetscFree(svd->data);
534: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
535: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
536: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
537: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
538: return(0);
539: }
543: PETSC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
544: {
546: SVD_CYCLIC *cyclic;
549: PetscNewLog(svd,&cyclic);
550: svd->data = (void*)cyclic;
551: svd->ops->solve = SVDSolve_Cyclic;
552: svd->ops->setup = SVDSetUp_Cyclic;
553: svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
554: svd->ops->destroy = SVDDestroy_Cyclic;
555: svd->ops->reset = SVDReset_Cyclic;
556: svd->ops->view = SVDView_Cyclic;
557: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
558: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
559: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
560: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
561: return(0);
562: }