Actual source code: cyclic.c
slepc-3.15.1 2021-05-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc singular value solver: "cyclic"
13: Method: Uses a Hermitian eigensolver for H(A) = [ 0 A ; A^T 0 ]
14: */
16: #include <slepc/private/svdimpl.h>
17: #include <slepc/private/epsimpl.h>
18: #include "cyclic.h"
20: static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
21: {
22: PetscErrorCode ierr;
23: SVD svd;
24: SVD_CYCLIC *cyclic;
25: const PetscScalar *px;
26: PetscScalar *py;
27: PetscInt m;
30: MatShellGetContext(B,(void**)&svd);
31: cyclic = (SVD_CYCLIC*)svd->data;
32: MatGetLocalSize(svd->A,&m,NULL);
33: VecGetArrayRead(x,&px);
34: VecGetArray(y,&py);
35: VecPlaceArray(cyclic->x1,px);
36: VecPlaceArray(cyclic->x2,px+m);
37: VecPlaceArray(cyclic->y1,py);
38: VecPlaceArray(cyclic->y2,py+m);
39: MatMult(svd->A,cyclic->x2,cyclic->y1);
40: MatMult(svd->AT,cyclic->x1,cyclic->y2);
41: VecResetArray(cyclic->x1);
42: VecResetArray(cyclic->x2);
43: VecResetArray(cyclic->y1);
44: VecResetArray(cyclic->y2);
45: VecRestoreArrayRead(x,&px);
46: VecRestoreArray(y,&py);
47: return(0);
48: }
50: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
51: {
55: VecSet(diag,0.0);
56: return(0);
57: }
59: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
60: {
61: PetscErrorCode ierr;
62: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
63: PetscInt M,N,m,n,i,isl,Istart,Iend;
64: const PetscScalar *isa;
65: PetscScalar *va;
66: PetscBool trackall,issinv;
67: Vec v;
68: Mat Zm,Zn;
69: ST st;
70: #if defined(PETSC_HAVE_CUDA)
71: PetscBool cuda;
72: #endif
75: SVDCheckStandard(svd);
76: MatGetSize(svd->A,&M,&N);
77: MatGetLocalSize(svd->A,&m,&n);
78: MatDestroy(&cyclic->mat);
79: if (cyclic->explicitmatrix) {
80: if (!svd->expltrans) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
81: MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
82: MatSetSizes(Zm,m,m,M,M);
83: MatSetFromOptions(Zm);
84: MatSetUp(Zm);
85: MatGetOwnershipRange(Zm,&Istart,&Iend);
86: for (i=Istart;i<Iend;i++) {
87: MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
88: }
89: MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
91: MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
92: MatSetSizes(Zn,n,n,N,N);
93: MatSetFromOptions(Zn);
94: MatSetUp(Zn);
95: MatGetOwnershipRange(Zn,&Istart,&Iend);
96: for (i=Istart;i<Iend;i++) {
97: MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
98: }
99: MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
101: MatCreateTile(1.0,Zm,1.0,svd->A,1.0,svd->AT,1.0,Zn,&cyclic->mat);
102: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->mat);
103: MatDestroy(&Zm);
104: MatDestroy(&Zn);
105: } else {
106: VecDestroy(&cyclic->x1);
107: VecDestroy(&cyclic->x2);
108: VecDestroy(&cyclic->y1);
109: VecDestroy(&cyclic->y2);
110: MatCreateVecsEmpty(svd->A,&cyclic->x2,&cyclic->x1);
111: MatCreateVecsEmpty(svd->A,&cyclic->y2,&cyclic->y1);
112: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->x1);
113: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->x2);
114: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->y1);
115: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->y2);
116: MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,svd,&cyclic->mat);
117: MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic);
118: #if defined(PETSC_HAVE_CUDA)
119: PetscObjectTypeCompareAny((PetscObject)(svd->A?svd->A:svd->AT),&cuda,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
120: if (cuda) {
121: MatShellSetVecType(cyclic->mat,VECCUDA);
122: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))MatMult_Cyclic_CUDA);
123: } else
124: #endif
125: {
126: MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))MatMult_Cyclic);
127: }
128: }
129: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->mat);
131: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
132: EPSSetOperators(cyclic->eps,cyclic->mat,NULL);
133: EPSSetProblemType(cyclic->eps,EPS_HEP);
134: if (!cyclic->usereps) {
135: if (svd->which == SVD_LARGEST) {
136: EPSGetST(cyclic->eps,&st);
137: PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv);
138: if (issinv) {
139: EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE);
140: } else {
141: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
142: }
143: } else {
144: EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
145: EPSSetTarget(cyclic->eps,0.0);
146: }
147: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
148: EPSSetTolerances(cyclic->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it);
149: switch (svd->conv) {
150: case SVD_CONV_ABS:
151: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_ABS);break;
152: case SVD_CONV_REL:
153: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_REL);break;
154: case SVD_CONV_MAXIT:
155: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
156: case SVD_CONV_USER:
157: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
158: }
159: }
160: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
161: /* Transfer the trackall option from svd to eps */
162: SVDGetTrackAll(svd,&trackall);
163: EPSSetTrackAll(cyclic->eps,trackall);
164: /* Transfer the initial subspace from svd to eps */
165: if (svd->nini<0 || svd->ninil<0) {
166: for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
167: MatCreateVecs(cyclic->mat,&v,NULL);
168: VecGetArray(v,&va);
169: if (i<-svd->ninil) {
170: VecGetSize(svd->ISL[i],&isl);
171: if (isl!=m) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
172: VecGetArrayRead(svd->ISL[i],&isa);
173: PetscArraycpy(va,isa,m);
174: VecRestoreArrayRead(svd->IS[i],&isa);
175: } else {
176: PetscArrayzero(&va,m);
177: }
178: if (i<-svd->nini) {
179: VecGetSize(svd->IS[i],&isl);
180: if (isl!=n) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
181: VecGetArrayRead(svd->IS[i],&isa);
182: PetscArraycpy(va+m,isa,n);
183: VecRestoreArrayRead(svd->IS[i],&isa);
184: } else {
185: PetscArrayzero(va+m,n);
186: }
187: VecRestoreArray(v,&va);
188: VecDestroy(&svd->IS[i]);
189: svd->IS[i] = v;
190: }
191: svd->nini = PetscMin(svd->nini,svd->ninil);
192: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
193: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
194: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
195: }
196: EPSSetUp(cyclic->eps);
197: EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
198: svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
199: EPSGetTolerances(cyclic->eps,NULL,&svd->max_it);
200: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
202: svd->leftbasis = PETSC_TRUE;
203: SVDAllocateSolution(svd,0);
204: return(0);
205: }
207: PetscErrorCode SVDSolve_Cyclic(SVD svd)
208: {
210: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
211: PetscInt i,j,nconv;
212: PetscScalar sigma;
215: EPSSolve(cyclic->eps);
216: EPSGetConverged(cyclic->eps,&nconv);
217: EPSGetIterationNumber(cyclic->eps,&svd->its);
218: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
219: for (i=0,j=0;i<nconv;i++) {
220: EPSGetEigenvalue(cyclic->eps,i,&sigma,NULL);
221: if (PetscRealPart(sigma) > 0.0) {
222: svd->sigma[j] = PetscRealPart(sigma);
223: j++;
224: }
225: }
226: svd->nconv = j;
227: return(0);
228: }
230: PetscErrorCode SVDComputeVectors_Cyclic(SVD svd)
231: {
232: PetscErrorCode ierr;
233: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
234: PetscInt i,j,M,N,m,n,nconv;
235: PetscScalar sigma;
236: const PetscScalar *px;
237: Vec x,x1,x2;
240: EPSGetConverged(cyclic->eps,&nconv);
241: MatCreateVecs(cyclic->mat,&x,NULL);
242: MatGetSize(svd->A,&M,&N);
243: MatGetLocalSize(svd->A,&m,&n);
244: MatCreateVecsEmpty(svd->A,&x2,&x1);
245: for (i=0,j=0;i<nconv;i++) {
246: EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
247: if (PetscRealPart(sigma) > 0.0) {
248: VecGetArrayRead(x,&px);
249: VecPlaceArray(x1,px);
250: VecPlaceArray(x2,px+m);
251: BVInsertVec(svd->U,j,x1);
252: BVScaleColumn(svd->U,j,PETSC_SQRT2);
253: BVInsertVec(svd->V,j,x2);
254: BVScaleColumn(svd->V,j,PETSC_SQRT2);
255: VecResetArray(x1);
256: VecResetArray(x2);
257: VecRestoreArrayRead(x,&px);
258: j++;
259: }
260: }
261: VecDestroy(&x);
262: VecDestroy(&x1);
263: VecDestroy(&x2);
264: return(0);
265: }
267: static PetscErrorCode EPSMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
268: {
269: PetscInt i,j;
270: SVD svd = (SVD)ctx;
271: PetscScalar er,ei;
275: nconv = 0;
276: for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
277: er = eigr[i]; ei = eigi[i];
278: STBackTransform(eps->st,1,&er,&ei);
279: if (PetscRealPart(er) > 0.0) {
280: svd->sigma[j] = PetscRealPart(er);
281: svd->errest[j] = errest[i];
282: if (errest[i] && errest[i] < svd->tol) nconv++;
283: j++;
284: }
285: }
286: nest = j;
287: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
288: return(0);
289: }
291: PetscErrorCode SVDSetFromOptions_Cyclic(PetscOptionItems *PetscOptionsObject,SVD svd)
292: {
294: PetscBool set,val;
295: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
296: ST st;
299: PetscOptionsHead(PetscOptionsObject,"SVD Cyclic Options");
301: PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
302: if (set) { SVDCyclicSetExplicitMatrix(svd,val); }
304: PetscOptionsTail();
306: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
307: if (!cyclic->explicitmatrix && !cyclic->usereps) {
308: /* use as default an ST with shell matrix and Jacobi */
309: EPSGetST(cyclic->eps,&st);
310: STSetMatMode(st,ST_MATMODE_SHELL);
311: }
312: EPSSetFromOptions(cyclic->eps);
313: return(0);
314: }
316: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
317: {
318: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
321: if (cyclic->explicitmatrix != explicitmatrix) {
322: cyclic->explicitmatrix = explicitmatrix;
323: svd->state = SVD_STATE_INITIAL;
324: }
325: return(0);
326: }
328: /*@
329: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
330: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
332: Logically Collective on svd
334: Input Parameters:
335: + svd - singular value solver
336: - explicit - boolean flag indicating if H(A) is built explicitly
338: Options Database Key:
339: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
341: Level: advanced
343: .seealso: SVDCyclicGetExplicitMatrix()
344: @*/
345: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
346: {
352: PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
353: return(0);
354: }
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: }
365: /*@
366: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly.
368: Not Collective
370: Input Parameter:
371: . svd - singular value solver
373: Output Parameter:
374: . explicit - the mode flag
376: Level: advanced
378: .seealso: SVDCyclicSetExplicitMatrix()
379: @*/
380: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
381: {
387: PetscUseMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
388: return(0);
389: }
391: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
392: {
393: PetscErrorCode ierr;
394: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
397: PetscObjectReference((PetscObject)eps);
398: EPSDestroy(&cyclic->eps);
399: cyclic->eps = eps;
400: cyclic->usereps = PETSC_TRUE;
401: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
402: svd->state = SVD_STATE_INITIAL;
403: return(0);
404: }
406: /*@
407: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
408: singular value solver.
410: Collective on svd
412: Input Parameters:
413: + svd - singular value solver
414: - eps - the eigensolver object
416: Level: advanced
418: .seealso: SVDCyclicGetEPS()
419: @*/
420: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
421: {
428: PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
429: return(0);
430: }
432: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
433: {
435: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
438: if (!cyclic->eps) {
439: EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
440: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
441: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
442: EPSAppendOptionsPrefix(cyclic->eps,"svd_cyclic_");
443: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
444: PetscObjectSetOptions((PetscObject)cyclic->eps,((PetscObject)svd)->options);
445: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
446: EPSMonitorSet(cyclic->eps,EPSMonitor_Cyclic,svd,NULL);
447: }
448: *eps = cyclic->eps;
449: return(0);
450: }
452: /*@
453: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
454: to the singular value solver.
456: Not Collective
458: Input Parameter:
459: . svd - singular value solver
461: Output Parameter:
462: . eps - the eigensolver object
464: Level: advanced
466: .seealso: SVDCyclicSetEPS()
467: @*/
468: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
469: {
475: PetscUseMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
476: return(0);
477: }
479: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
480: {
482: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
483: PetscBool isascii;
486: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
487: if (isascii) {
488: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
489: PetscViewerASCIIPrintf(viewer," %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
490: PetscViewerASCIIPushTab(viewer);
491: EPSView(cyclic->eps,viewer);
492: PetscViewerASCIIPopTab(viewer);
493: }
494: return(0);
495: }
497: PetscErrorCode SVDReset_Cyclic(SVD svd)
498: {
500: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
503: EPSReset(cyclic->eps);
504: MatDestroy(&cyclic->mat);
505: VecDestroy(&cyclic->x1);
506: VecDestroy(&cyclic->x2);
507: VecDestroy(&cyclic->y1);
508: VecDestroy(&cyclic->y2);
509: return(0);
510: }
512: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
513: {
515: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
518: EPSDestroy(&cyclic->eps);
519: PetscFree(svd->data);
520: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
521: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
522: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
523: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
524: return(0);
525: }
527: SLEPC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
528: {
530: SVD_CYCLIC *cyclic;
533: PetscNewLog(svd,&cyclic);
534: svd->data = (void*)cyclic;
535: svd->ops->solve = SVDSolve_Cyclic;
536: svd->ops->setup = SVDSetUp_Cyclic;
537: svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
538: svd->ops->destroy = SVDDestroy_Cyclic;
539: svd->ops->reset = SVDReset_Cyclic;
540: svd->ops->view = SVDView_Cyclic;
541: svd->ops->computevectors = SVDComputeVectors_Cyclic;
542: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
543: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
544: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
545: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
546: return(0);
547: }