Actual source code: cross.c

slepc-3.19.0 2023-03-31
Report Typos and Errors
  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:     PetscCall(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: 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(MatSetUp(Omega));
191:     PetscCall(MatDiagonalSet(Omega,svd->omega,INSERT_VALUES));
192:     PetscCall(EPSSetOperators(cross->eps,cross->C,Omega));
193:     PetscCall(EPSSetProblemType(cross->eps,EPS_GHIEP));
194:     PetscCall(MatDestroy(&Omega));
195:   } else {
196:     PetscCall(EPSSetOperators(cross->eps,cross->C,NULL));
197:     PetscCall(EPSSetProblemType(cross->eps,EPS_HEP));
198:   }
199:   if (!cross->usereps) {
200:     PetscCall(EPSGetST(cross->eps,&st));
201:     PetscCall(PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv));
202:     PetscCall(PetscObjectTypeCompare((PetscObject)cross->eps,EPSKRYLOVSCHUR,&isks));
203:     if (svd->isgeneralized && svd->which==SVD_SMALLEST) {
204:       if (cross->explicitmatrix && isks && !issinv) {  /* default to shift-and-invert */
205:         PetscCall(STSetType(st,STSINVERT));
206:         PetscCall(EPSSetTarget(cross->eps,0.0));
207:         which = EPS_TARGET_REAL;
208:       } else which = issinv?EPS_TARGET_REAL:EPS_SMALLEST_REAL;
209:     } else {
210:       if (issinv) which = EPS_TARGET_MAGNITUDE;
211:       else if (svd->ishyperbolic) which = svd->which==SVD_LARGEST?EPS_LARGEST_MAGNITUDE:EPS_SMALLEST_MAGNITUDE;
212:       else which = svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL;
213:     }
214:     PetscCall(EPSSetWhichEigenpairs(cross->eps,which));
215:     PetscCall(EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd));
216:     PetscCall(EPSSetTolerances(cross->eps,svd->tol==(PetscReal)PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it));
217:     switch (svd->conv) {
218:     case SVD_CONV_ABS:
219:       PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS));break;
220:     case SVD_CONV_REL:
221:       PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_REL));break;
222:     case SVD_CONV_NORM:
223:       if (svd->isgeneralized) {
224:         if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
225:         if (!svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
226:         PetscCall(EPSSetConvergenceTestFunction(cross->eps,EPSConv_Cross,svd,NULL));
227:       } else {
228:         PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_NORM));break;
229:       }
230:       break;
231:     case SVD_CONV_MAXIT:
232:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
233:     case SVD_CONV_USER:
234:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
235:     }
236:   }
237:   SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
238:   /* Transfer the trackall option from svd to eps */
239:   PetscCall(SVDGetTrackAll(svd,&trackall));
240:   PetscCall(EPSSetTrackAll(cross->eps,trackall));
241:   /* Transfer the initial space from svd to eps */
242:   if (svd->nini<0) {
243:     PetscCall(EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS));
244:     PetscCall(SlepcBasisDestroy_Private(&svd->nini,&svd->IS));
245:   }
246:   PetscCall(EPSSetUp(cross->eps));
247:   PetscCall(EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd));
248:   PetscCall(EPSGetTolerances(cross->eps,NULL,&svd->max_it));
249:   if (svd->tol==(PetscReal)PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

251:   svd->leftbasis = PETSC_FALSE;
252:   PetscCall(SVDAllocateSolution(svd,0));
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

256: PetscErrorCode SVDSolve_Cross(SVD svd)
257: {
258:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
259:   PetscInt       i;
260:   PetscScalar    lambda;
261:   PetscReal      sigma;

263:   PetscFunctionBegin;
264:   PetscCall(EPSSolve(cross->eps));
265:   PetscCall(EPSGetConverged(cross->eps,&svd->nconv));
266:   PetscCall(EPSGetIterationNumber(cross->eps,&svd->its));
267:   PetscCall(EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason));
268:   for (i=0;i<svd->nconv;i++) {
269:     PetscCall(EPSGetEigenvalue(cross->eps,i,&lambda,NULL));
270:     sigma = PetscRealPart(lambda);
271:     if (svd->ishyperbolic) svd->sigma[i] = PetscSqrtReal(PetscAbsReal(sigma));
272:     else {
273:       PetscCheck(sigma>-10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_FP,"Negative eigenvalue computed by EPS: %g",(double)sigma);
274:       if (sigma<0.0) {
275:         PetscCall(PetscInfo(svd,"Negative eigenvalue computed by EPS: %g, resetting to 0\n",(double)sigma));
276:         sigma = 0.0;
277:       }
278:       svd->sigma[i] = PetscSqrtReal(sigma);
279:     }
280:   }
281:   PetscFunctionReturn(PETSC_SUCCESS);
282: }

284: PetscErrorCode SVDComputeVectors_Cross(SVD svd)
285: {
286:   SVD_CROSS         *cross = (SVD_CROSS*)svd->data;
287:   PetscInt          i,mloc,ploc;
288:   Vec               u,v,x,uv,w,omega2=NULL;
289:   Mat               Omega;
290:   PetscScalar       *dst,alpha,lambda,*varray;
291:   const PetscScalar *src;
292:   PetscReal         nrm;

294:   PetscFunctionBegin;
295:   if (svd->isgeneralized) {
296:     PetscCall(MatCreateVecs(svd->A,NULL,&u));
297:     PetscCall(VecGetLocalSize(u,&mloc));
298:     PetscCall(MatCreateVecs(svd->B,NULL,&v));
299:     PetscCall(VecGetLocalSize(v,&ploc));
300:     for (i=0;i<svd->nconv;i++) {
301:       PetscCall(BVGetColumn(svd->V,i,&x));
302:       PetscCall(EPSGetEigenpair(cross->eps,i,&lambda,NULL,x,NULL));
303:       PetscCall(MatMult(svd->A,x,u));     /* u_i*c_i/alpha = A*x_i */
304:       PetscCall(VecNormalize(u,NULL));
305:       PetscCall(MatMult(svd->B,x,v));     /* v_i*s_i/alpha = B*x_i */
306:       PetscCall(VecNormalize(v,&nrm));    /* ||v||_2 = s_i/alpha   */
307:       alpha = 1.0/(PetscSqrtReal(1.0+PetscRealPart(lambda))*nrm);    /* alpha=s_i/||v||_2 */
308:       PetscCall(VecScale(x,alpha));
309:       PetscCall(BVRestoreColumn(svd->V,i,&x));
310:       /* copy [u;v] to U[i] */
311:       PetscCall(BVGetColumn(svd->U,i,&uv));
312:       PetscCall(VecGetArrayWrite(uv,&dst));
313:       PetscCall(VecGetArrayRead(u,&src));
314:       PetscCall(PetscArraycpy(dst,src,mloc));
315:       PetscCall(VecRestoreArrayRead(u,&src));
316:       PetscCall(VecGetArrayRead(v,&src));
317:       PetscCall(PetscArraycpy(dst+mloc,src,ploc));
318:       PetscCall(VecRestoreArrayRead(v,&src));
319:       PetscCall(VecRestoreArrayWrite(uv,&dst));
320:       PetscCall(BVRestoreColumn(svd->U,i,&uv));
321:     }
322:     PetscCall(VecDestroy(&v));
323:     PetscCall(VecDestroy(&u));
324:   } else if (svd->ishyperbolic && svd->swapped) {  /* was solved as GHIEP, set u=Omega*u and normalize */
325:     PetscCall(EPSGetOperators(cross->eps,NULL,&Omega));
326:     PetscCall(MatCreateVecs(Omega,&w,NULL));
327:     PetscCall(VecCreateSeq(PETSC_COMM_SELF,svd->ncv,&omega2));
328:     PetscCall(VecGetArrayWrite(omega2,&varray));
329:     for (i=0;i<svd->nconv;i++) {
330:       PetscCall(BVGetColumn(svd->V,i,&v));
331:       PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
332:       PetscCall(MatMult(Omega,v,w));
333:       PetscCall(VecDot(v,w,&alpha));
334:       svd->sign[i] = PetscSign(PetscRealPart(alpha));
335:       varray[i] = svd->sign[i];
336:       alpha = 1.0/PetscSqrtScalar(PetscAbsScalar(alpha));
337:       PetscCall(VecScale(w,alpha));
338:       PetscCall(VecCopy(w,v));
339:       PetscCall(BVRestoreColumn(svd->V,i,&v));
340:     }
341:     PetscCall(BVSetSignature(svd->V,omega2));
342:     PetscCall(VecRestoreArrayWrite(omega2,&varray));
343:     PetscCall(VecDestroy(&omega2));
344:     PetscCall(VecDestroy(&w));
345:     PetscCall(SVDComputeVectors_Left(svd));
346:   } else {
347:     for (i=0;i<svd->nconv;i++) {
348:       PetscCall(BVGetColumn(svd->V,i,&v));
349:       PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
350:       PetscCall(BVRestoreColumn(svd->V,i,&v));
351:     }
352:     PetscCall(SVDComputeVectors_Left(svd));
353:   }
354:   PetscFunctionReturn(PETSC_SUCCESS);
355: }

357: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
358: {
359:   PetscInt       i;
360:   SVD            svd = (SVD)ctx;
361:   PetscScalar    er,ei;
362:   ST             st;

364:   PetscFunctionBegin;
365:   PetscCall(EPSGetST(eps,&st));
366:   for (i=0;i<PetscMin(nest,svd->ncv);i++) {
367:     er = eigr[i]; ei = eigi[i];
368:     PetscCall(STBackTransform(st,1,&er,&ei));
369:     svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
370:     svd->errest[i] = errest[i];
371:   }
372:   PetscCall(SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest));
373:   PetscFunctionReturn(PETSC_SUCCESS);
374: }

376: PetscErrorCode SVDSetFromOptions_Cross(SVD svd,PetscOptionItems *PetscOptionsObject)
377: {
378:   PetscBool      set,val;
379:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
380:   ST             st;

382:   PetscFunctionBegin;
383:   PetscOptionsHeadBegin(PetscOptionsObject,"SVD Cross Options");

385:     PetscCall(PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set));
386:     if (set) PetscCall(SVDCrossSetExplicitMatrix(svd,val));

388:   PetscOptionsHeadEnd();

390:   if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
391:   if (!cross->explicitmatrix && !cross->usereps) {
392:     /* use as default an ST with shell matrix and Jacobi */
393:     PetscCall(EPSGetST(cross->eps,&st));
394:     PetscCall(STSetMatMode(st,ST_MATMODE_SHELL));
395:   }
396:   PetscCall(EPSSetFromOptions(cross->eps));
397:   PetscFunctionReturn(PETSC_SUCCESS);
398: }

400: static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
401: {
402:   SVD_CROSS *cross = (SVD_CROSS*)svd->data;

404:   PetscFunctionBegin;
405:   if (cross->explicitmatrix != explicitmatrix) {
406:     cross->explicitmatrix = explicitmatrix;
407:     svd->state = SVD_STATE_INITIAL;
408:   }
409:   PetscFunctionReturn(PETSC_SUCCESS);
410: }

412: /*@
413:    SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
414:    be computed explicitly.

416:    Logically Collective

418:    Input Parameters:
419: +  svd         - singular value solver
420: -  explicitmat - boolean flag indicating if A^T*A is built explicitly

422:    Options Database Key:
423: .  -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag

425:    Level: advanced

427: .seealso: SVDCrossGetExplicitMatrix()
428: @*/
429: PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmat)
430: {
431:   PetscFunctionBegin;
434:   PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmat)
439: {
440:   SVD_CROSS *cross = (SVD_CROSS*)svd->data;

442:   PetscFunctionBegin;
443:   *explicitmat = cross->explicitmatrix;
444:   PetscFunctionReturn(PETSC_SUCCESS);
445: }

447: /*@
448:    SVDCrossGetExplicitMatrix - Returns the flag indicating if A^T*A is built explicitly.

450:    Not Collective

452:    Input Parameter:
453: .  svd  - singular value solver

455:    Output Parameter:
456: .  explicitmat - the mode flag

458:    Level: advanced

460: .seealso: SVDCrossSetExplicitMatrix()
461: @*/
462: PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
463: {
464:   PetscFunctionBegin;
467:   PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
468:   PetscFunctionReturn(PETSC_SUCCESS);
469: }

471: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
472: {
473:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

475:   PetscFunctionBegin;
476:   PetscCall(PetscObjectReference((PetscObject)eps));
477:   PetscCall(EPSDestroy(&cross->eps));
478:   cross->eps     = eps;
479:   cross->usereps = PETSC_TRUE;
480:   svd->state     = SVD_STATE_INITIAL;
481:   PetscFunctionReturn(PETSC_SUCCESS);
482: }

484: /*@
485:    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
486:    singular value solver.

488:    Collective

490:    Input Parameters:
491: +  svd - singular value solver
492: -  eps - the eigensolver object

494:    Level: advanced

496: .seealso: SVDCrossGetEPS()
497: @*/
498: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
499: {
500:   PetscFunctionBegin;
503:   PetscCheckSameComm(svd,1,eps,2);
504:   PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }

508: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
509: {
510:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

512:   PetscFunctionBegin;
513:   if (!cross->eps) {
514:     PetscCall(EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps));
515:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1));
516:     PetscCall(EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix));
517:     PetscCall(EPSAppendOptionsPrefix(cross->eps,"svd_cross_"));
518:     PetscCall(PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options));
519:     PetscCall(EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL));
520:     PetscCall(EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL));
521:   }
522:   *eps = cross->eps;
523:   PetscFunctionReturn(PETSC_SUCCESS);
524: }

526: /*@
527:    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
528:    to the singular value solver.

530:    Collective

532:    Input Parameter:
533: .  svd - singular value solver

535:    Output Parameter:
536: .  eps - the eigensolver object

538:    Level: advanced

540: .seealso: SVDCrossSetEPS()
541: @*/
542: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
543: {
544:   PetscFunctionBegin;
547:   PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
548:   PetscFunctionReturn(PETSC_SUCCESS);
549: }

551: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
552: {
553:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
554:   PetscBool      isascii;

556:   PetscFunctionBegin;
557:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
558:   if (isascii) {
559:     if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
560:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %s matrix\n",cross->explicitmatrix?"explicit":"implicit"));
561:     PetscCall(PetscViewerASCIIPushTab(viewer));
562:     PetscCall(EPSView(cross->eps,viewer));
563:     PetscCall(PetscViewerASCIIPopTab(viewer));
564:   }
565:   PetscFunctionReturn(PETSC_SUCCESS);
566: }

568: PetscErrorCode SVDReset_Cross(SVD svd)
569: {
570:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

572:   PetscFunctionBegin;
573:   PetscCall(EPSReset(cross->eps));
574:   PetscCall(MatDestroy(&cross->C));
575:   PetscCall(MatDestroy(&cross->D));
576:   PetscFunctionReturn(PETSC_SUCCESS);
577: }

579: PetscErrorCode SVDDestroy_Cross(SVD svd)
580: {
581:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

583:   PetscFunctionBegin;
584:   PetscCall(EPSDestroy(&cross->eps));
585:   PetscCall(PetscFree(svd->data));
586:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL));
587:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL));
588:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL));
589:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL));
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
594: {
595:   SVD_CROSS      *cross;

597:   PetscFunctionBegin;
598:   PetscCall(PetscNew(&cross));
599:   svd->data = (void*)cross;

601:   svd->ops->solve          = SVDSolve_Cross;
602:   svd->ops->solveg         = SVDSolve_Cross;
603:   svd->ops->solveh         = SVDSolve_Cross;
604:   svd->ops->setup          = SVDSetUp_Cross;
605:   svd->ops->setfromoptions = SVDSetFromOptions_Cross;
606:   svd->ops->destroy        = SVDDestroy_Cross;
607:   svd->ops->reset          = SVDReset_Cross;
608:   svd->ops->view           = SVDView_Cross;
609:   svd->ops->computevectors = SVDComputeVectors_Cross;
610:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross));
611:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross));
612:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross));
613:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross));
614:   PetscFunctionReturn(PETSC_SUCCESS);
615: }