Actual source code: cross.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: "cross"
13: Method: Uses a Hermitian eigensolver for A^T*A
14: */
16: #include <slepc/private/svdimpl.h>
17: #include <slepc/private/epsimpl.h>
19: typedef struct {
20: PetscBool explicitmatrix;
21: EPS eps;
22: PetscBool usereps;
23: Mat mat;
24: Vec w,diag;
25: } SVD_CROSS;
27: static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
28: {
30: SVD svd;
31: SVD_CROSS *cross;
34: MatShellGetContext(B,(void**)&svd);
35: cross = (SVD_CROSS*)svd->data;
36: MatMult(svd->A,x,cross->w);
37: MatMult(svd->AT,cross->w,y);
38: return(0);
39: }
41: static PetscErrorCode MatCreateVecs_Cross(Mat B,Vec *right,Vec *left)
42: {
44: SVD svd;
47: MatShellGetContext(B,(void**)&svd);
48: if (right) {
49: MatCreateVecs(svd->A,right,NULL);
50: if (left) { VecDuplicate(*right,left); }
51: } else {
52: MatCreateVecs(svd->A,left,NULL);
53: }
54: return(0);
55: }
57: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
58: {
59: PetscErrorCode ierr;
60: SVD svd;
61: SVD_CROSS *cross;
62: PetscMPIInt len;
63: PetscInt N,n,i,j,start,end,ncols;
64: PetscScalar *work1,*work2,*diag;
65: const PetscInt *cols;
66: const PetscScalar *vals;
69: MatShellGetContext(B,(void**)&svd);
70: cross = (SVD_CROSS*)svd->data;
71: if (!cross->diag) {
72: /* compute diagonal from rows and store in cross->diag */
73: VecDuplicate(d,&cross->diag);
74: MatGetSize(svd->A,NULL,&N);
75: MatGetLocalSize(svd->A,NULL,&n);
76: PetscCalloc2(N,&work1,N,&work2);
77: if (svd->swapped) {
78: MatGetOwnershipRange(svd->AT,&start,&end);
79: for (i=start;i<end;i++) {
80: MatGetRow(svd->AT,i,&ncols,NULL,&vals);
81: for (j=0;j<ncols;j++)
82: work1[i] += vals[j]*vals[j];
83: MatRestoreRow(svd->AT,i,&ncols,NULL,&vals);
84: }
85: } else {
86: MatGetOwnershipRange(svd->A,&start,&end);
87: for (i=start;i<end;i++) {
88: MatGetRow(svd->A,i,&ncols,&cols,&vals);
89: for (j=0;j<ncols;j++)
90: work1[cols[j]] += vals[j]*vals[j];
91: MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
92: }
93: }
94: PetscMPIIntCast(N,&len);
95: MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)svd));
96: VecGetOwnershipRange(cross->diag,&start,&end);
97: VecGetArray(cross->diag,&diag);
98: for (i=start;i<end;i++) diag[i-start] = work2[i];
99: VecRestoreArray(cross->diag,&diag);
100: PetscFree2(work1,work2);
101: }
102: VecCopy(cross->diag,d);
103: return(0);
104: }
106: PetscErrorCode SVDSetUp_Cross(SVD svd)
107: {
109: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
110: PetscInt n;
111: PetscBool trackall;
114: SVDCheckStandard(svd);
115: MatDestroy(&cross->mat);
116: VecDestroy(&cross->w);
117: if (cross->explicitmatrix) {
118: if (svd->expltrans) { /* explicit transpose */
119: MatProductCreate(svd->AT,svd->A,NULL,&cross->mat);
120: MatProductSetType(cross->mat,MATPRODUCT_AB);
121: } else { /* implicit transpose */
122: #if defined(PETSC_USE_COMPLEX)
123: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Must use explicit transpose with complex scalars");
124: #else
125: if (!svd->swapped) {
126: MatProductCreate(svd->A,svd->A,NULL,&cross->mat);
127: MatProductSetType(cross->mat,MATPRODUCT_AtB);
128: } else {
129: MatProductCreate(svd->AT,svd->AT,NULL,&cross->mat);
130: MatProductSetType(cross->mat,MATPRODUCT_ABt);
131: }
132: #endif
133: }
134: MatProductSetFromOptions(cross->mat);
135: MatProductSymbolic(cross->mat);
136: MatProductNumeric(cross->mat);
137: } else {
138: MatGetLocalSize(svd->A,NULL,&n);
139: MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
140: MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
141: MatShellSetOperation(cross->mat,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_Cross);
142: MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
143: MatCreateVecs(svd->A,NULL,&cross->w);
144: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->w);
145: }
146: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->mat);
148: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
149: EPSSetOperators(cross->eps,cross->mat,NULL);
150: EPSSetProblemType(cross->eps,EPS_HEP);
151: if (!cross->usereps) {
152: EPSSetWhichEigenpairs(cross->eps,svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL);
153: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd);
154: EPSSetTolerances(cross->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it);
155: switch (svd->conv) {
156: case SVD_CONV_ABS:
157: EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS);break;
158: case SVD_CONV_REL:
159: EPSSetConvergenceTest(cross->eps,EPS_CONV_REL);break;
160: case SVD_CONV_MAXIT:
161: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
162: case SVD_CONV_USER:
163: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
164: }
165: }
166: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
167: /* Transfer the trackall option from svd to eps */
168: SVDGetTrackAll(svd,&trackall);
169: EPSSetTrackAll(cross->eps,trackall);
170: /* Transfer the initial space from svd to eps */
171: if (svd->nini<0) {
172: EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
173: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
174: }
175: EPSSetUp(cross->eps);
176: EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
177: EPSGetTolerances(cross->eps,NULL,&svd->max_it);
178: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
180: svd->leftbasis = PETSC_FALSE;
181: SVDAllocateSolution(svd,0);
182: return(0);
183: }
185: PetscErrorCode SVDSolve_Cross(SVD svd)
186: {
188: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
189: PetscInt i;
190: PetscScalar lambda;
191: PetscReal sigma;
194: EPSSolve(cross->eps);
195: EPSGetConverged(cross->eps,&svd->nconv);
196: EPSGetIterationNumber(cross->eps,&svd->its);
197: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
198: for (i=0;i<svd->nconv;i++) {
199: EPSGetEigenvalue(cross->eps,i,&lambda,NULL);
200: sigma = PetscRealPart(lambda);
201: if (sigma<-10*PETSC_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS: %g",sigma);
202: if (sigma<0.0) {
203: PetscInfo1(svd,"Negative eigenvalue computed by EPS: %g, resetting to 0\n",sigma);
204: sigma = 0.0;
205: }
206: svd->sigma[i] = PetscSqrtReal(sigma);
207: }
208: return(0);
209: }
211: PetscErrorCode SVDComputeVectors_Cross(SVD svd)
212: {
214: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
215: PetscInt i;
216: Vec v;
219: for (i=0;i<svd->nconv;i++) {
220: BVGetColumn(svd->V,i,&v);
221: EPSGetEigenvector(cross->eps,i,v,NULL);
222: BVRestoreColumn(svd->V,i,&v);
223: }
224: SVDComputeVectors_Left(svd);
225: return(0);
226: }
228: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
229: {
230: PetscInt i;
231: SVD svd = (SVD)ctx;
232: PetscScalar er,ei;
236: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
237: er = eigr[i]; ei = eigi[i];
238: STBackTransform(eps->st,1,&er,&ei);
239: svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
240: svd->errest[i] = errest[i];
241: }
242: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
243: return(0);
244: }
246: PetscErrorCode SVDSetFromOptions_Cross(PetscOptionItems *PetscOptionsObject,SVD svd)
247: {
249: PetscBool set,val;
250: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
251: ST st;
254: PetscOptionsHead(PetscOptionsObject,"SVD Cross Options");
256: PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set);
257: if (set) { SVDCrossSetExplicitMatrix(svd,val); }
259: PetscOptionsTail();
261: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
262: if (!cross->explicitmatrix && !cross->usereps) {
263: /* use as default an ST with shell matrix and Jacobi */
264: EPSGetST(cross->eps,&st);
265: STSetMatMode(st,ST_MATMODE_SHELL);
266: }
267: EPSSetFromOptions(cross->eps);
268: return(0);
269: }
271: static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
272: {
273: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
276: if (cross->explicitmatrix != explicitmatrix) {
277: cross->explicitmatrix = explicitmatrix;
278: svd->state = SVD_STATE_INITIAL;
279: }
280: return(0);
281: }
283: /*@
284: SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
285: be computed explicitly.
287: Logically Collective on svd
289: Input Parameters:
290: + svd - singular value solver
291: - explicit - boolean flag indicating if A^T*A is built explicitly
293: Options Database Key:
294: . -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag
296: Level: advanced
298: .seealso: SVDCrossGetExplicitMatrix()
299: @*/
300: PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
301: {
307: PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
308: return(0);
309: }
311: static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmatrix)
312: {
313: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
316: *explicitmatrix = cross->explicitmatrix;
317: return(0);
318: }
320: /*@
321: SVDCrossGetExplicitMatrix - Returns the flag indicating if A^T*A is built explicitly.
323: Not Collective
325: Input Parameter:
326: . svd - singular value solver
328: Output Parameter:
329: . explicit - the mode flag
331: Level: advanced
333: .seealso: SVDCrossSetExplicitMatrix()
334: @*/
335: PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
336: {
342: PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
343: return(0);
344: }
346: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
347: {
349: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
352: PetscObjectReference((PetscObject)eps);
353: EPSDestroy(&cross->eps);
354: cross->eps = eps;
355: cross->usereps = PETSC_TRUE;
356: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
357: svd->state = SVD_STATE_INITIAL;
358: return(0);
359: }
361: /*@
362: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
363: singular value solver.
365: Collective on svd
367: Input Parameters:
368: + svd - singular value solver
369: - eps - the eigensolver object
371: Level: advanced
373: .seealso: SVDCrossGetEPS()
374: @*/
375: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
376: {
383: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
384: return(0);
385: }
387: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
388: {
389: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
393: if (!cross->eps) {
394: EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
395: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
396: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
397: EPSAppendOptionsPrefix(cross->eps,"svd_cross_");
398: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
399: PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options);
400: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
401: EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL);
402: }
403: *eps = cross->eps;
404: return(0);
405: }
407: /*@
408: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
409: to the singular value solver.
411: Not Collective
413: Input Parameter:
414: . svd - singular value solver
416: Output Parameter:
417: . eps - the eigensolver object
419: Level: advanced
421: .seealso: SVDCrossSetEPS()
422: @*/
423: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
424: {
430: PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
431: return(0);
432: }
434: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
435: {
437: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
438: PetscBool isascii;
441: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
442: if (isascii) {
443: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
444: PetscViewerASCIIPrintf(viewer," %s matrix\n",cross->explicitmatrix?"explicit":"implicit");
445: PetscViewerASCIIPushTab(viewer);
446: EPSView(cross->eps,viewer);
447: PetscViewerASCIIPopTab(viewer);
448: }
449: return(0);
450: }
452: PetscErrorCode SVDReset_Cross(SVD svd)
453: {
455: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
458: EPSReset(cross->eps);
459: MatDestroy(&cross->mat);
460: VecDestroy(&cross->w);
461: VecDestroy(&cross->diag);
462: return(0);
463: }
465: PetscErrorCode SVDDestroy_Cross(SVD svd)
466: {
468: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
471: EPSDestroy(&cross->eps);
472: PetscFree(svd->data);
473: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
474: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
475: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL);
476: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL);
477: return(0);
478: }
480: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
481: {
483: SVD_CROSS *cross;
486: PetscNewLog(svd,&cross);
487: svd->data = (void*)cross;
489: svd->ops->solve = SVDSolve_Cross;
490: svd->ops->setup = SVDSetUp_Cross;
491: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
492: svd->ops->destroy = SVDDestroy_Cross;
493: svd->ops->reset = SVDReset_Cross;
494: svd->ops->view = SVDView_Cross;
495: svd->ops->computevectors = SVDComputeVectors_Cross;
496: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
497: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
498: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross);
499: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross);
500: return(0);
501: }