Actual source code: cross.c
slepc-main 2024-11-15
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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>
18: typedef struct {
19: PetscBool explicitmatrix;
20: EPS eps;
21: PetscBool usereps;
22: Mat C,D;
23: } SVD_CROSS;
25: typedef struct {
26: Mat A,AT;
27: Vec w,diag,omega;
28: PetscBool swapped;
29: } SVD_CROSS_SHELL;
31: static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
32: {
33: SVD_CROSS_SHELL *ctx;
35: PetscFunctionBegin;
36: PetscCall(MatShellGetContext(B,&ctx));
37: PetscCall(MatMult(ctx->A,x,ctx->w));
38: if (ctx->omega && !ctx->swapped) PetscCall(VecPointwiseMult(ctx->w,ctx->w,ctx->omega));
39: PetscCall(MatMult(ctx->AT,ctx->w,y));
40: PetscFunctionReturn(PETSC_SUCCESS);
41: }
43: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
44: {
45: SVD_CROSS_SHELL *ctx;
46: PetscMPIInt len;
47: PetscInt N,n,i,j,start,end,ncols;
48: PetscScalar *work1,*work2,*diag;
49: const PetscInt *cols;
50: const PetscScalar *vals;
52: PetscFunctionBegin;
53: PetscCall(MatShellGetContext(B,&ctx));
54: if (!ctx->diag) {
55: /* compute diagonal from rows and store in ctx->diag */
56: PetscCall(VecDuplicate(d,&ctx->diag));
57: PetscCall(MatGetSize(ctx->A,NULL,&N));
58: PetscCall(MatGetLocalSize(ctx->A,NULL,&n));
59: PetscCall(PetscCalloc2(N,&work1,N,&work2));
60: if (ctx->swapped) {
61: PetscCall(MatGetOwnershipRange(ctx->AT,&start,&end));
62: for (i=start;i<end;i++) {
63: PetscCall(MatGetRow(ctx->AT,i,&ncols,NULL,&vals));
64: for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
65: PetscCall(MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals));
66: }
67: } else {
68: PetscCall(MatGetOwnershipRange(ctx->A,&start,&end));
69: for (i=start;i<end;i++) {
70: PetscCall(MatGetRow(ctx->A,i,&ncols,&cols,&vals));
71: for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
72: PetscCall(MatRestoreRow(ctx->A,i,&ncols,&cols,&vals));
73: }
74: }
75: PetscCall(PetscMPIIntCast(N,&len));
76: PetscCallMPI(MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B)));
77: PetscCall(VecGetOwnershipRange(ctx->diag,&start,&end));
78: PetscCall(VecGetArrayWrite(ctx->diag,&diag));
79: for (i=start;i<end;i++) diag[i-start] = work2[i];
80: PetscCall(VecRestoreArrayWrite(ctx->diag,&diag));
81: PetscCall(PetscFree2(work1,work2));
82: }
83: PetscCall(VecCopy(ctx->diag,d));
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: static PetscErrorCode MatDestroy_Cross(Mat B)
88: {
89: SVD_CROSS_SHELL *ctx;
91: PetscFunctionBegin;
92: PetscCall(MatShellGetContext(B,&ctx));
93: PetscCall(VecDestroy(&ctx->w));
94: PetscCall(VecDestroy(&ctx->diag));
95: PetscCall(PetscFree(ctx));
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: static PetscErrorCode SVDCrossGetProductMat(SVD svd,Mat A,Mat AT,Mat *C)
100: {
101: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
102: SVD_CROSS_SHELL *ctx;
103: PetscInt n;
104: VecType vtype;
105: Mat B;
107: PetscFunctionBegin;
108: if (cross->explicitmatrix) {
109: if (!svd->ishyperbolic || svd->swapped) B = (!svd->expltrans && svd->swapped)? AT: A;
110: else { /* duplicate A and scale by signature */
111: PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
112: PetscCall(MatDiagonalScale(B,svd->omega,NULL));
113: }
114: if (svd->expltrans) { /* explicit transpose */
115: PetscCall(MatProductCreate(AT,B,NULL,C));
116: PetscCall(MatProductSetType(*C,MATPRODUCT_AB));
117: } else { /* implicit transpose */
118: #if defined(PETSC_USE_COMPLEX)
119: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Must use explicit transpose with complex scalars");
120: #else
121: if (!svd->swapped) {
122: PetscCall(MatProductCreate(A,B,NULL,C));
123: PetscCall(MatProductSetType(*C,MATPRODUCT_AtB));
124: } else {
125: PetscCall(MatProductCreate(B,AT,NULL,C));
126: PetscCall(MatProductSetType(*C,MATPRODUCT_ABt));
127: }
128: #endif
129: }
130: PetscCall(MatProductSetFromOptions(*C));
131: PetscCall(MatProductSymbolic(*C));
132: PetscCall(MatProductNumeric(*C));
133: if (svd->ishyperbolic && !svd->swapped) PetscCall(MatDestroy(&B));
134: } else {
135: PetscCall(PetscNew(&ctx));
136: ctx->A = A;
137: ctx->AT = AT;
138: ctx->omega = svd->omega;
139: ctx->swapped = svd->swapped;
140: PetscCall(MatCreateVecs(A,NULL,&ctx->w));
141: PetscCall(MatGetLocalSize(A,NULL,&n));
142: PetscCall(MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,(void*)ctx,C));
143: PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cross));
144: if (!svd->ishyperbolic || svd->swapped) PetscCall(MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross));
145: PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cross));
146: PetscCall(MatGetVecType(A,&vtype));
147: PetscCall(MatSetVecType(*C,vtype));
148: }
149: PetscFunctionReturn(PETSC_SUCCESS);
150: }
152: /* Convergence test relative to the norm of R (used in GSVD only) */
153: static PetscErrorCode EPSConv_Cross(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
154: {
155: SVD svd = (SVD)ctx;
157: PetscFunctionBegin;
158: *errest = res/PetscMax(svd->nrma,svd->nrmb);
159: PetscFunctionReturn(PETSC_SUCCESS);
160: }
162: static PetscErrorCode SVDSetUp_Cross(SVD svd)
163: {
164: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
165: ST st;
166: PetscBool trackall,issinv,isks;
167: EPSProblemType ptype;
168: EPSWhich which;
169: Mat Omega;
170: MatType Atype;
171: PetscInt n,N;
173: PetscFunctionBegin;
174: if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
175: PetscCall(MatDestroy(&cross->C));
176: PetscCall(MatDestroy(&cross->D));
177: PetscCall(SVDCrossGetProductMat(svd,svd->A,svd->AT,&cross->C));
178: if (svd->isgeneralized) {
179: PetscCall(SVDCrossGetProductMat(svd,svd->B,svd->BT,&cross->D));
180: PetscCall(EPSSetOperators(cross->eps,cross->C,cross->D));
181: PetscCall(EPSGetProblemType(cross->eps,&ptype));
182: if (!ptype) PetscCall(EPSSetProblemType(cross->eps,EPS_GHEP));
183: } else if (svd->ishyperbolic && svd->swapped) {
184: PetscCall(MatGetType(svd->OP,&Atype));
185: PetscCall(MatGetSize(svd->A,NULL,&N));
186: PetscCall(MatGetLocalSize(svd->A,NULL,&n));
187: PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Omega));
188: PetscCall(MatSetSizes(Omega,n,n,N,N));
189: PetscCall(MatSetType(Omega,Atype));
190: PetscCall(MatDiagonalSet(Omega,svd->omega,INSERT_VALUES));
191: PetscCall(EPSSetOperators(cross->eps,cross->C,Omega));
192: PetscCall(EPSSetProblemType(cross->eps,EPS_GHIEP));
193: PetscCall(MatDestroy(&Omega));
194: } else {
195: PetscCall(EPSSetOperators(cross->eps,cross->C,NULL));
196: PetscCall(EPSSetProblemType(cross->eps,EPS_HEP));
197: }
198: if (!cross->usereps) {
199: PetscCall(EPSGetST(cross->eps,&st));
200: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv));
201: PetscCall(PetscObjectTypeCompare((PetscObject)cross->eps,EPSKRYLOVSCHUR,&isks));
202: if (svd->isgeneralized && svd->which==SVD_SMALLEST) {
203: if (cross->explicitmatrix && isks && !issinv) { /* default to shift-and-invert */
204: PetscCall(STSetType(st,STSINVERT));
205: PetscCall(EPSSetTarget(cross->eps,0.0));
206: which = EPS_TARGET_REAL;
207: } else which = issinv?EPS_TARGET_REAL:EPS_SMALLEST_REAL;
208: } else {
209: if (issinv) which = EPS_TARGET_MAGNITUDE;
210: else if (svd->ishyperbolic) which = svd->which==SVD_LARGEST?EPS_LARGEST_MAGNITUDE:EPS_SMALLEST_MAGNITUDE;
211: else which = svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL;
212: }
213: PetscCall(EPSSetWhichEigenpairs(cross->eps,which));
214: PetscCall(EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd));
215: PetscCall(EPSSetTolerances(cross->eps,svd->tol==(PetscReal)PETSC_DETERMINE?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it));
216: switch (svd->conv) {
217: case SVD_CONV_ABS:
218: PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS));break;
219: case SVD_CONV_REL:
220: PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_REL));break;
221: case SVD_CONV_NORM:
222: if (svd->isgeneralized) {
223: if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
224: if (!svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
225: PetscCall(EPSSetConvergenceTestFunction(cross->eps,EPSConv_Cross,svd,NULL));
226: } else {
227: PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_NORM));break;
228: }
229: break;
230: case SVD_CONV_MAXIT:
231: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
232: case SVD_CONV_USER:
233: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
234: }
235: }
236: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
237: /* Transfer the trackall option from svd to eps */
238: PetscCall(SVDGetTrackAll(svd,&trackall));
239: PetscCall(EPSSetTrackAll(cross->eps,trackall));
240: /* Transfer the initial space from svd to eps */
241: if (svd->nini<0) {
242: PetscCall(EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS));
243: PetscCall(SlepcBasisDestroy_Private(&svd->nini,&svd->IS));
244: }
245: PetscCall(EPSSetUp(cross->eps));
246: PetscCall(EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd));
247: PetscCall(EPSGetTolerances(cross->eps,NULL,&svd->max_it));
248: if (svd->tol==(PetscReal)PETSC_DETERMINE) svd->tol = SLEPC_DEFAULT_TOL;
250: svd->leftbasis = PETSC_FALSE;
251: PetscCall(SVDAllocateSolution(svd,0));
252: PetscFunctionReturn(PETSC_SUCCESS);
253: }
255: static PetscErrorCode SVDSolve_Cross(SVD svd)
256: {
257: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
258: PetscInt i;
259: PetscScalar lambda;
260: PetscReal sigma;
262: PetscFunctionBegin;
263: PetscCall(EPSSolve(cross->eps));
264: PetscCall(EPSGetConverged(cross->eps,&svd->nconv));
265: PetscCall(EPSGetIterationNumber(cross->eps,&svd->its));
266: PetscCall(EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason));
267: for (i=0;i<svd->nconv;i++) {
268: PetscCall(EPSGetEigenvalue(cross->eps,i,&lambda,NULL));
269: sigma = PetscRealPart(lambda);
270: if (svd->ishyperbolic) svd->sigma[i] = PetscSqrtReal(PetscAbsReal(sigma));
271: else {
272: PetscCheck(sigma>-10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_FP,"Negative eigenvalue computed by EPS: %g",(double)sigma);
273: if (sigma<0.0) {
274: PetscCall(PetscInfo(svd,"Negative eigenvalue computed by EPS: %g, resetting to 0\n",(double)sigma));
275: sigma = 0.0;
276: }
277: svd->sigma[i] = PetscSqrtReal(sigma);
278: }
279: }
280: PetscFunctionReturn(PETSC_SUCCESS);
281: }
283: static PetscErrorCode SVDComputeVectors_Cross(SVD svd)
284: {
285: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
286: PetscInt i,mloc,ploc;
287: Vec u,v,x,uv,w,omega2=NULL;
288: Mat Omega;
289: PetscScalar *dst,alpha,lambda,*varray;
290: const PetscScalar *src;
291: PetscReal nrm;
293: PetscFunctionBegin;
294: if (svd->isgeneralized) {
295: PetscCall(MatCreateVecs(svd->A,NULL,&u));
296: PetscCall(VecGetLocalSize(u,&mloc));
297: PetscCall(MatCreateVecs(svd->B,NULL,&v));
298: PetscCall(VecGetLocalSize(v,&ploc));
299: for (i=0;i<svd->nconv;i++) {
300: PetscCall(BVGetColumn(svd->V,i,&x));
301: PetscCall(EPSGetEigenpair(cross->eps,i,&lambda,NULL,x,NULL));
302: PetscCall(MatMult(svd->A,x,u)); /* u_i*c_i/alpha = A*x_i */
303: PetscCall(VecNormalize(u,NULL));
304: PetscCall(MatMult(svd->B,x,v)); /* v_i*s_i/alpha = B*x_i */
305: PetscCall(VecNormalize(v,&nrm)); /* ||v||_2 = s_i/alpha */
306: alpha = 1.0/(PetscSqrtReal(1.0+PetscRealPart(lambda))*nrm); /* alpha=s_i/||v||_2 */
307: PetscCall(VecScale(x,alpha));
308: PetscCall(BVRestoreColumn(svd->V,i,&x));
309: /* copy [u;v] to U[i] */
310: PetscCall(BVGetColumn(svd->U,i,&uv));
311: PetscCall(VecGetArrayWrite(uv,&dst));
312: PetscCall(VecGetArrayRead(u,&src));
313: PetscCall(PetscArraycpy(dst,src,mloc));
314: PetscCall(VecRestoreArrayRead(u,&src));
315: PetscCall(VecGetArrayRead(v,&src));
316: PetscCall(PetscArraycpy(dst+mloc,src,ploc));
317: PetscCall(VecRestoreArrayRead(v,&src));
318: PetscCall(VecRestoreArrayWrite(uv,&dst));
319: PetscCall(BVRestoreColumn(svd->U,i,&uv));
320: }
321: PetscCall(VecDestroy(&v));
322: PetscCall(VecDestroy(&u));
323: } else if (svd->ishyperbolic && svd->swapped) { /* was solved as GHIEP, set u=Omega*u and normalize */
324: PetscCall(EPSGetOperators(cross->eps,NULL,&Omega));
325: PetscCall(MatCreateVecs(Omega,&w,NULL));
326: PetscCall(VecCreateSeq(PETSC_COMM_SELF,svd->ncv,&omega2));
327: PetscCall(VecGetArrayWrite(omega2,&varray));
328: for (i=0;i<svd->nconv;i++) {
329: PetscCall(BVGetColumn(svd->V,i,&v));
330: PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
331: PetscCall(MatMult(Omega,v,w));
332: PetscCall(VecDot(v,w,&alpha));
333: svd->sign[i] = PetscSign(PetscRealPart(alpha));
334: varray[i] = svd->sign[i];
335: alpha = 1.0/PetscSqrtScalar(PetscAbsScalar(alpha));
336: PetscCall(VecScale(w,alpha));
337: PetscCall(VecCopy(w,v));
338: PetscCall(BVRestoreColumn(svd->V,i,&v));
339: }
340: PetscCall(BVSetSignature(svd->V,omega2));
341: PetscCall(VecRestoreArrayWrite(omega2,&varray));
342: PetscCall(VecDestroy(&omega2));
343: PetscCall(VecDestroy(&w));
344: PetscCall(SVDComputeVectors_Left(svd));
345: } else {
346: for (i=0;i<svd->nconv;i++) {
347: PetscCall(BVGetColumn(svd->V,i,&v));
348: PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
349: PetscCall(BVRestoreColumn(svd->V,i,&v));
350: }
351: PetscCall(SVDComputeVectors_Left(svd));
352: }
353: PetscFunctionReturn(PETSC_SUCCESS);
354: }
356: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
357: {
358: PetscInt i;
359: SVD svd = (SVD)ctx;
360: PetscScalar er,ei;
361: ST st;
363: PetscFunctionBegin;
364: PetscCall(EPSGetST(eps,&st));
365: for (i=0;i<PetscMin(nest,svd->ncv);i++) {
366: er = eigr[i]; ei = eigi[i];
367: PetscCall(STBackTransform(st,1,&er,&ei));
368: svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
369: svd->errest[i] = errest[i];
370: }
371: PetscCall(SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest));
372: PetscFunctionReturn(PETSC_SUCCESS);
373: }
375: static PetscErrorCode SVDSetFromOptions_Cross(SVD svd,PetscOptionItems *PetscOptionsObject)
376: {
377: PetscBool set,val;
378: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
379: ST st;
381: PetscFunctionBegin;
382: PetscOptionsHeadBegin(PetscOptionsObject,"SVD Cross Options");
384: PetscCall(PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set));
385: if (set) PetscCall(SVDCrossSetExplicitMatrix(svd,val));
387: PetscOptionsHeadEnd();
389: if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
390: if (!cross->explicitmatrix && !cross->usereps) {
391: /* use as default an ST with shell matrix and Jacobi */
392: PetscCall(EPSGetST(cross->eps,&st));
393: PetscCall(STSetMatMode(st,ST_MATMODE_SHELL));
394: }
395: PetscCall(EPSSetFromOptions(cross->eps));
396: PetscFunctionReturn(PETSC_SUCCESS);
397: }
399: static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
400: {
401: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
403: PetscFunctionBegin;
404: if (cross->explicitmatrix != explicitmatrix) {
405: cross->explicitmatrix = explicitmatrix;
406: svd->state = SVD_STATE_INITIAL;
407: }
408: PetscFunctionReturn(PETSC_SUCCESS);
409: }
411: /*@
412: SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
413: be computed explicitly.
415: Logically Collective
417: Input Parameters:
418: + svd - singular value solver
419: - explicitmat - boolean flag indicating if A^T*A is built explicitly
421: Options Database Key:
422: . -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag
424: Level: advanced
426: .seealso: SVDCrossGetExplicitMatrix()
427: @*/
428: PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmat)
429: {
430: PetscFunctionBegin;
433: PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
434: PetscFunctionReturn(PETSC_SUCCESS);
435: }
437: static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmat)
438: {
439: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
441: PetscFunctionBegin;
442: *explicitmat = cross->explicitmatrix;
443: PetscFunctionReturn(PETSC_SUCCESS);
444: }
446: /*@
447: SVDCrossGetExplicitMatrix - Returns the flag indicating if A^T*A is built explicitly.
449: Not Collective
451: Input Parameter:
452: . svd - singular value solver
454: Output Parameter:
455: . explicitmat - the mode flag
457: Level: advanced
459: .seealso: SVDCrossSetExplicitMatrix()
460: @*/
461: PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
462: {
463: PetscFunctionBegin;
465: PetscAssertPointer(explicitmat,2);
466: PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
467: PetscFunctionReturn(PETSC_SUCCESS);
468: }
470: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
471: {
472: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
474: PetscFunctionBegin;
475: PetscCall(PetscObjectReference((PetscObject)eps));
476: PetscCall(EPSDestroy(&cross->eps));
477: cross->eps = eps;
478: cross->usereps = PETSC_TRUE;
479: svd->state = SVD_STATE_INITIAL;
480: PetscFunctionReturn(PETSC_SUCCESS);
481: }
483: /*@
484: SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
485: singular value solver.
487: Collective
489: Input Parameters:
490: + svd - singular value solver
491: - eps - the eigensolver object
493: Level: advanced
495: .seealso: SVDCrossGetEPS()
496: @*/
497: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
498: {
499: PetscFunctionBegin;
502: PetscCheckSameComm(svd,1,eps,2);
503: PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
508: {
509: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
511: PetscFunctionBegin;
512: if (!cross->eps) {
513: PetscCall(EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps));
514: PetscCall(PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1));
515: PetscCall(EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix));
516: PetscCall(EPSAppendOptionsPrefix(cross->eps,"svd_cross_"));
517: PetscCall(PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options));
518: PetscCall(EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL));
519: PetscCall(EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL));
520: }
521: *eps = cross->eps;
522: PetscFunctionReturn(PETSC_SUCCESS);
523: }
525: /*@
526: SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
527: to the singular value solver.
529: Collective
531: Input Parameter:
532: . svd - singular value solver
534: Output Parameter:
535: . eps - the eigensolver object
537: Level: advanced
539: .seealso: SVDCrossSetEPS()
540: @*/
541: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
542: {
543: PetscFunctionBegin;
545: PetscAssertPointer(eps,2);
546: PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
547: PetscFunctionReturn(PETSC_SUCCESS);
548: }
550: static PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
551: {
552: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
553: PetscBool isascii;
555: PetscFunctionBegin;
556: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
557: if (isascii) {
558: if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
559: PetscCall(PetscViewerASCIIPrintf(viewer," %s matrix\n",cross->explicitmatrix?"explicit":"implicit"));
560: PetscCall(PetscViewerASCIIPushTab(viewer));
561: PetscCall(EPSView(cross->eps,viewer));
562: PetscCall(PetscViewerASCIIPopTab(viewer));
563: }
564: PetscFunctionReturn(PETSC_SUCCESS);
565: }
567: static PetscErrorCode SVDReset_Cross(SVD svd)
568: {
569: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
571: PetscFunctionBegin;
572: PetscCall(EPSReset(cross->eps));
573: PetscCall(MatDestroy(&cross->C));
574: PetscCall(MatDestroy(&cross->D));
575: PetscFunctionReturn(PETSC_SUCCESS);
576: }
578: static PetscErrorCode SVDDestroy_Cross(SVD svd)
579: {
580: SVD_CROSS *cross = (SVD_CROSS*)svd->data;
582: PetscFunctionBegin;
583: PetscCall(EPSDestroy(&cross->eps));
584: PetscCall(PetscFree(svd->data));
585: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL));
586: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL));
587: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL));
588: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL));
589: PetscFunctionReturn(PETSC_SUCCESS);
590: }
592: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
593: {
594: SVD_CROSS *cross;
596: PetscFunctionBegin;
597: PetscCall(PetscNew(&cross));
598: svd->data = (void*)cross;
600: svd->ops->solve = SVDSolve_Cross;
601: svd->ops->solveg = SVDSolve_Cross;
602: svd->ops->solveh = SVDSolve_Cross;
603: svd->ops->setup = SVDSetUp_Cross;
604: svd->ops->setfromoptions = SVDSetFromOptions_Cross;
605: svd->ops->destroy = SVDDestroy_Cross;
606: svd->ops->reset = SVDReset_Cross;
607: svd->ops->view = SVDView_Cross;
608: svd->ops->computevectors = SVDComputeVectors_Cross;
609: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross));
610: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross));
611: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross));
612: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross));
613: PetscFunctionReturn(PETSC_SUCCESS);
614: }