Actual source code: cross.c
slepc-3.16.2 2022-02-01
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 C,D;
24: } SVD_CROSS;
26: typedef struct {
27: Mat A,AT;
28: Vec w,diag;
29: PetscBool swapped;
30: } SVD_CROSS_SHELL;
32: static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
33: {
34: PetscErrorCode ierr;
35: SVD_CROSS_SHELL *ctx;
38: MatShellGetContext(B,&ctx);
39: MatMult(ctx->A,x,ctx->w);
40: MatMult(ctx->AT,ctx->w,y);
41: return(0);
42: }
44: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
45: {
46: PetscErrorCode ierr;
47: SVD_CROSS_SHELL *ctx;
48: PetscMPIInt len;
49: PetscInt N,n,i,j,start,end,ncols;
50: PetscScalar *work1,*work2,*diag;
51: const PetscInt *cols;
52: const PetscScalar *vals;
55: MatShellGetContext(B,&ctx);
56: if (!ctx->diag) {
57: /* compute diagonal from rows and store in ctx->diag */
58: VecDuplicate(d,&ctx->diag);
59: MatGetSize(ctx->A,NULL,&N);
60: MatGetLocalSize(ctx->A,NULL,&n);
61: PetscCalloc2(N,&work1,N,&work2);
62: if (ctx->swapped) {
63: MatGetOwnershipRange(ctx->AT,&start,&end);
64: for (i=start;i<end;i++) {
65: MatGetRow(ctx->AT,i,&ncols,NULL,&vals);
66: for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
67: MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals);
68: }
69: } else {
70: MatGetOwnershipRange(ctx->A,&start,&end);
71: for (i=start;i<end;i++) {
72: MatGetRow(ctx->A,i,&ncols,&cols,&vals);
73: for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
74: MatRestoreRow(ctx->A,i,&ncols,&cols,&vals);
75: }
76: }
77: PetscMPIIntCast(N,&len);
78: MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B));
79: VecGetOwnershipRange(ctx->diag,&start,&end);
80: VecGetArrayWrite(ctx->diag,&diag);
81: for (i=start;i<end;i++) diag[i-start] = work2[i];
82: VecRestoreArrayWrite(ctx->diag,&diag);
83: PetscFree2(work1,work2);
84: }
85: VecCopy(ctx->diag,d);
86: return(0);
87: }
89: static PetscErrorCode MatDestroy_Cross(Mat B)
90: {
91: PetscErrorCode ierr;
92: SVD_CROSS_SHELL *ctx;
95: MatShellGetContext(B,&ctx);
96: VecDestroy(&ctx->w);
97: VecDestroy(&ctx->diag);
98: PetscFree(ctx);
99: return(0);
100: }
102: static PetscErrorCode SVDCrossGetProductMat(SVD svd,Mat A,Mat AT,Mat *C)
103: {
104: PetscErrorCode ierr;
105: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
106: SVD_CROSS_SHELL *ctx;
107: PetscInt n;
108: VecType vtype;
111: if (cross->explicitmatrix) {
112: if (svd->expltrans) { /* explicit transpose */
113: MatProductCreate(AT,A,NULL,C);
114: MatProductSetType(*C,MATPRODUCT_AB);
115: } else { /* implicit transpose */
116: #if defined(PETSC_USE_COMPLEX)
117: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Must use explicit transpose with complex scalars");
118: #else
119: if (!svd->swapped) {
120: MatProductCreate(A,A,NULL,C);
121: MatProductSetType(*C,MATPRODUCT_AtB);
122: } else {
123: MatProductCreate(AT,AT,NULL,C);
124: MatProductSetType(*C,MATPRODUCT_ABt);
125: }
126: #endif
127: }
128: MatProductSetFromOptions(*C);
129: MatProductSymbolic(*C);
130: MatProductNumeric(*C);
131: } else {
132: PetscNew(&ctx);
133: ctx->A = A;
134: ctx->AT = AT;
135: ctx->swapped = svd->swapped;
136: MatCreateVecs(A,NULL,&ctx->w);
137: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->w);
138: MatGetLocalSize(A,NULL,&n);
139: MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,(void*)ctx,C);
140: MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cross);
141: MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
142: MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cross);
143: MatGetVecType(A,&vtype);
144: MatSetVecType(*C,vtype);
145: }
146: PetscLogObjectParent((PetscObject)svd,(PetscObject)*C);
147: return(0);
148: }
150: PetscErrorCode SVDSetUp_Cross(SVD svd)
151: {
153: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
154: ST st;
155: PetscBool trackall,issinv;
158: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
159: MatDestroy(&cross->C);
160: MatDestroy(&cross->D);
161: SVDCrossGetProductMat(svd,svd->A,svd->AT,&cross->C);
162: if (svd->isgeneralized) {
163: SVDCrossGetProductMat(svd,svd->B,svd->BT,&cross->D);
164: EPSSetOperators(cross->eps,cross->C,cross->D);
165: EPSSetProblemType(cross->eps,EPS_GHEP);
166: } else {
167: EPSSetOperators(cross->eps,cross->C,NULL);
168: EPSSetProblemType(cross->eps,EPS_HEP);
169: }
170: if (!cross->usereps) {
171: EPSGetST(cross->eps,&st);
172: if (svd->isgeneralized && svd->which==SVD_SMALLEST) {
173: STSetType(st,STSINVERT);
174: EPSSetTarget(cross->eps,0.0);
175: EPSSetWhichEigenpairs(cross->eps,EPS_TARGET_REAL);
176: } else {
177: PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv);
178: if (issinv) {
179: EPSSetWhichEigenpairs(cross->eps,EPS_TARGET_MAGNITUDE);
180: } else {
181: EPSSetWhichEigenpairs(cross->eps,svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL);
182: }
183: }
184: EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd);
185: EPSSetTolerances(cross->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it);
186: switch (svd->conv) {
187: case SVD_CONV_ABS:
188: EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS);break;
189: case SVD_CONV_REL:
190: EPSSetConvergenceTest(cross->eps,EPS_CONV_REL);break;
191: case SVD_CONV_NORM:
192: EPSSetConvergenceTest(cross->eps,EPS_CONV_NORM);break;
193: case SVD_CONV_MAXIT:
194: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
195: case SVD_CONV_USER:
196: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
197: }
198: }
199: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
200: /* Transfer the trackall option from svd to eps */
201: SVDGetTrackAll(svd,&trackall);
202: EPSSetTrackAll(cross->eps,trackall);
203: /* Transfer the initial space from svd to eps */
204: if (svd->nini<0) {
205: EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
206: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
207: }
208: EPSSetUp(cross->eps);
209: EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
210: EPSGetTolerances(cross->eps,NULL,&svd->max_it);
211: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
213: svd->leftbasis = PETSC_FALSE;
214: SVDAllocateSolution(svd,0);
215: return(0);
216: }
218: PetscErrorCode SVDSolve_Cross(SVD svd)
219: {
221: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
222: PetscInt i;
223: PetscScalar lambda;
224: PetscReal sigma;
227: EPSSolve(cross->eps);
228: EPSGetConverged(cross->eps,&svd->nconv);
229: EPSGetIterationNumber(cross->eps,&svd->its);
230: EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
231: for (i=0;i<svd->nconv;i++) {
232: EPSGetEigenvalue(cross->eps,i,&lambda,NULL);
233: sigma = PetscRealPart(lambda);
234: if (sigma<-10*PETSC_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS: %g",sigma);
235: if (sigma<0.0) {
236: PetscInfo1(svd,"Negative eigenvalue computed by EPS: %g, resetting to 0\n",sigma);
237: sigma = 0.0;
238: }
239: svd->sigma[i] = PetscSqrtReal(sigma);
240: }
241: return(0);
242: }
244: PetscErrorCode SVDComputeVectors_Cross(SVD svd)
245: {
246: PetscErrorCode ierr;
247: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
248: PetscInt i,mloc,ploc;
249: Vec u,v,x,uv;
250: PetscScalar *dst,alpha,lambda;
251: const PetscScalar *src;
252: PetscReal nrm;
255: if (svd->isgeneralized) {
256: MatCreateVecs(svd->A,NULL,&u);
257: VecGetLocalSize(u,&mloc);
258: MatCreateVecs(svd->B,NULL,&v);
259: VecGetLocalSize(v,&ploc);
260: for (i=0;i<svd->nconv;i++) {
261: BVGetColumn(svd->V,i,&x);
262: EPSGetEigenpair(cross->eps,i,&lambda,NULL,x,NULL);
263: MatMult(svd->A,x,u); /* u_i*c_i/alpha = A*x_i */
264: VecNormalize(u,NULL);
265: MatMult(svd->B,x,v); /* v_i*s_i/alpha = B*x_i */
266: VecNormalize(v,&nrm); /* ||v||_2 = s_i/alpha */
267: alpha = 1.0/(PetscSqrtReal(1.0+PetscRealPart(lambda))*nrm); /* alpha=s_i/||v||_2 */
268: VecScale(x,alpha);
269: BVRestoreColumn(svd->V,i,&x);
270: /* copy [u;v] to U[i] */
271: BVGetColumn(svd->U,i,&uv);
272: VecGetArrayWrite(uv,&dst);
273: VecGetArrayRead(u,&src);
274: PetscArraycpy(dst,src,mloc);
275: VecRestoreArrayRead(u,&src);
276: VecGetArrayRead(v,&src);
277: PetscArraycpy(dst+mloc,src,ploc);
278: VecRestoreArrayRead(v,&src);
279: VecRestoreArrayWrite(uv,&dst);
280: BVRestoreColumn(svd->U,i,&uv);
281: }
282: VecDestroy(&v);
283: VecDestroy(&u);
284: } else {
285: for (i=0;i<svd->nconv;i++) {
286: BVGetColumn(svd->V,i,&v);
287: EPSGetEigenvector(cross->eps,i,v,NULL);
288: BVRestoreColumn(svd->V,i,&v);
289: }
290: SVDComputeVectors_Left(svd);
291: }
292: return(0);
293: }
295: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
296: {
297: PetscInt i;
298: SVD svd = (SVD)ctx;
299: PetscScalar er,ei;
303: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
304: er = eigr[i]; ei = eigi[i];
305: STBackTransform(eps->st,1,&er,&ei);
306: svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
307: svd->errest[i] = errest[i];
308: }
309: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
310: return(0);
311: }
313: PetscErrorCode SVDSetFromOptions_Cross(PetscOptionItems *PetscOptionsObject,SVD svd)
314: {
316: PetscBool set,val;
317: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
318: ST st;
321: PetscOptionsHead(PetscOptionsObject,"SVD Cross Options");
323: PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set);
324: if (set) { SVDCrossSetExplicitMatrix(svd,val); }
326: PetscOptionsTail();
328: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
329: if (!cross->explicitmatrix && !cross->usereps) {
330: /* use as default an ST with shell matrix and Jacobi */
331: EPSGetST(cross->eps,&st);
332: STSetMatMode(st,ST_MATMODE_SHELL);
333: }
334: EPSSetFromOptions(cross->eps);
335: return(0);
336: }
338: static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
339: {
340: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
343: if (cross->explicitmatrix != explicitmatrix) {
344: cross->explicitmatrix = explicitmatrix;
345: svd->state = SVD_STATE_INITIAL;
346: }
347: return(0);
348: }
350: /*@
351: SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
352: be computed explicitly.
354: Logically Collective on svd
356: Input Parameters:
357: + svd - singular value solver
358: - explicit - boolean flag indicating if A^T*A is built explicitly
360: Options Database Key:
361: . -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag
363: Level: advanced
365: .seealso: SVDCrossGetExplicitMatrix()
366: @*/
367: PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
368: {
374: PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
375: return(0);
376: }
378: static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmatrix)
379: {
380: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
383: *explicitmatrix = cross->explicitmatrix;
384: return(0);
385: }
387: /*@
388: SVDCrossGetExplicitMatrix - Returns the flag indicating if A^T*A is built explicitly.
390: Not Collective
392: Input Parameter:
393: . svd - singular value solver
395: Output Parameter:
396: . explicit - the mode flag
398: Level: advanced
400: .seealso: SVDCrossSetExplicitMatrix()
401: @*/
402: PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
403: {
409: PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
410: return(0);
411: }
413: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
414: {
416: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
419: PetscObjectReference((PetscObject)eps);
420: EPSDestroy(&cross->eps);
421: cross->eps = eps;
422: cross->usereps = PETSC_TRUE;
423: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
424: svd->state = SVD_STATE_INITIAL;
425: return(0);
426: }
428: /*@
429: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
430: singular value solver.
432: Collective on svd
434: Input Parameters:
435: + svd - singular value solver
436: - eps - the eigensolver object
438: Level: advanced
440: .seealso: SVDCrossGetEPS()
441: @*/
442: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
443: {
450: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
451: return(0);
452: }
454: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
455: {
456: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
460: if (!cross->eps) {
461: EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
462: PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
463: EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
464: EPSAppendOptionsPrefix(cross->eps,"svd_cross_");
465: PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
466: PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options);
467: EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
468: EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL);
469: }
470: *eps = cross->eps;
471: return(0);
472: }
474: /*@
475: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
476: to the singular value solver.
478: Not Collective
480: Input Parameter:
481: . svd - singular value solver
483: Output Parameter:
484: . eps - the eigensolver object
486: Level: advanced
488: .seealso: SVDCrossSetEPS()
489: @*/
490: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
491: {
497: PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
498: return(0);
499: }
501: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
502: {
504: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
505: PetscBool isascii;
508: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
509: if (isascii) {
510: if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
511: PetscViewerASCIIPrintf(viewer," %s matrix\n",cross->explicitmatrix?"explicit":"implicit");
512: PetscViewerASCIIPushTab(viewer);
513: EPSView(cross->eps,viewer);
514: PetscViewerASCIIPopTab(viewer);
515: }
516: return(0);
517: }
519: PetscErrorCode SVDReset_Cross(SVD svd)
520: {
522: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
525: EPSReset(cross->eps);
526: MatDestroy(&cross->C);
527: MatDestroy(&cross->D);
528: return(0);
529: }
531: PetscErrorCode SVDDestroy_Cross(SVD svd)
532: {
534: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
537: EPSDestroy(&cross->eps);
538: PetscFree(svd->data);
539: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
540: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
541: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL);
542: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL);
543: return(0);
544: }
546: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
547: {
549: SVD_CROSS *cross;
552: PetscNewLog(svd,&cross);
553: svd->data = (void*)cross;
555: svd->ops->solve = SVDSolve_Cross;
556: svd->ops->solveg = SVDSolve_Cross;
557: svd->ops->setup = SVDSetUp_Cross;
558: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
559: svd->ops->destroy = SVDDestroy_Cross;
560: svd->ops->reset = SVDReset_Cross;
561: svd->ops->view = SVDView_Cross;
562: svd->ops->computevectors = SVDComputeVectors_Cross;
563: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
564: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
565: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross);
566: PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross);
567: return(0);
568: }