Actual source code: cross.c

slepc-3.13.2 2020-05-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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:   SVDMatMult(svd,PETSC_FALSE,x,cross->w);
 37:   SVDMatMult(svd,PETSC_TRUE,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:     SVDMatCreateVecs(svd,right,NULL);
 50:     if (left) { VecDuplicate(*right,left); }
 51:   } else {
 52:     SVDMatCreateVecs(svd,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:     SVDMatGetSize(svd,NULL,&N);
 75:     SVDMatGetLocalSize(svd,NULL,&n);
 76:     PetscCalloc2(N,&work1,N,&work2);
 77:     if (svd->AT) {
 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:     MPI_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:   if (!cross->mat) {
115:     if (cross->explicitmatrix) {
116:       if (svd->A && svd->AT) {  /* explicit transpose */
117:         MatProductCreate(svd->AT,svd->A,NULL,&cross->mat);
118:         MatProductSetType(cross->mat,MATPRODUCT_AB);
119:       } else {  /* implicit transpose */
120: #if defined(PETSC_USE_COMPLEX)
121:         SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Must use explicit transpose with complex scalars");
122: #else
123:         if (svd->A) {
124:           MatProductCreate(svd->A,svd->A,NULL,&cross->mat);
125:           MatProductSetType(cross->mat,MATPRODUCT_AtB);
126:         } else {
127:           MatProductCreate(svd->AT,svd->AT,NULL,&cross->mat);
128:           MatProductSetType(cross->mat,MATPRODUCT_ABt);
129:         }
130: #endif
131:       }
132:       MatProductSetFromOptions(cross->mat);
133:       MatProductSymbolic(cross->mat);
134:       MatProductNumeric(cross->mat);
135:     } else {
136:       SVDMatGetLocalSize(svd,NULL,&n);
137:       MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
138:       MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
139:       MatShellSetOperation(cross->mat,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_Cross);
140:       MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
141:       SVDMatCreateVecs(svd,NULL,&cross->w);
142:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->w);
143:     }
144:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->mat);
145:   }

147:   if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
148:   EPSSetOperators(cross->eps,cross->mat,NULL);
149:   EPSSetProblemType(cross->eps,EPS_HEP);
150:   if (!cross->usereps) {
151:     EPSSetWhichEigenpairs(cross->eps,svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL);
152:     EPSSetDimensions(cross->eps,svd->nsv,svd->ncv?svd->ncv:PETSC_DEFAULT,svd->mpd?svd->mpd:PETSC_DEFAULT);
153:     EPSSetTolerances(cross->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it?svd->max_it:PETSC_DEFAULT);
154:     switch (svd->conv) {
155:     case SVD_CONV_ABS:
156:       EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS);break;
157:     case SVD_CONV_REL:
158:       EPSSetConvergenceTest(cross->eps,EPS_CONV_REL);break;
159:     case SVD_CONV_USER:
160:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
161:     }
162:   }
163:   if (svd->stop!=SVD_STOP_BASIC) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined stopping test not supported in this solver");
164:   /* Transfer the trackall option from svd to eps */
165:   SVDGetTrackAll(svd,&trackall);
166:   EPSSetTrackAll(cross->eps,trackall);
167:   /* Transfer the initial space from svd to eps */
168:   if (svd->nini<0) {
169:     EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
170:     SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
171:   }
172:   EPSSetUp(cross->eps);
173:   EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
174:   EPSGetTolerances(cross->eps,NULL,&svd->max_it);
175:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

177:   svd->leftbasis = PETSC_FALSE;
178:   SVDAllocateSolution(svd,0);
179:   return(0);
180: }

182: PetscErrorCode SVDSolve_Cross(SVD svd)
183: {
185:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
186:   PetscInt       i;
187:   PetscScalar    lambda;
188:   PetscReal      sigma;
189:   Vec            v;

192:   EPSSolve(cross->eps);
193:   EPSGetConverged(cross->eps,&svd->nconv);
194:   EPSGetIterationNumber(cross->eps,&svd->its);
195:   EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
196:   for (i=0;i<svd->nconv;i++) {
197:     BVGetColumn(svd->V,i,&v);
198:     EPSGetEigenpair(cross->eps,i,&lambda,NULL,v,NULL);
199:     BVRestoreColumn(svd->V,i,&v);
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: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
212: {
213:   PetscInt       i;
214:   SVD            svd = (SVD)ctx;
215:   PetscScalar    er,ei;

219:   for (i=0;i<PetscMin(nest,svd->ncv);i++) {
220:     er = eigr[i]; ei = eigi[i];
221:     STBackTransform(eps->st,1,&er,&ei);
222:     svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
223:     svd->errest[i] = errest[i];
224:   }
225:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
226:   return(0);
227: }

229: PetscErrorCode SVDSetFromOptions_Cross(PetscOptionItems *PetscOptionsObject,SVD svd)
230: {
232:   PetscBool      set,val;
233:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
234:   ST             st;

237:   PetscOptionsHead(PetscOptionsObject,"SVD Cross Options");

239:     PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set);
240:     if (set) { SVDCrossSetExplicitMatrix(svd,val); }

242:   PetscOptionsTail();

244:   if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
245:   if (!cross->explicitmatrix && !cross->usereps) {
246:     /* use as default an ST with shell matrix and Jacobi */
247:     EPSGetST(cross->eps,&st);
248:     STSetMatMode(st,ST_MATMODE_SHELL);
249:   }
250:   EPSSetFromOptions(cross->eps);
251:   return(0);
252: }

254: static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
255: {
256:   SVD_CROSS *cross = (SVD_CROSS*)svd->data;

259:   if (cross->explicitmatrix != explicitmatrix) {
260:     cross->explicitmatrix = explicitmatrix;
261:     svd->state = SVD_STATE_INITIAL;
262:   }
263:   return(0);
264: }

266: /*@
267:    SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
268:    be computed explicitly.

270:    Logically Collective on svd

272:    Input Parameters:
273: +  svd      - singular value solver
274: -  explicit - boolean flag indicating if A^T*A is built explicitly

276:    Options Database Key:
277: .  -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag

279:    Level: advanced

281: .seealso: SVDCrossGetExplicitMatrix()
282: @*/
283: PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
284: {

290:   PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
291:   return(0);
292: }

294: static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmatrix)
295: {
296:   SVD_CROSS *cross = (SVD_CROSS*)svd->data;

299:   *explicitmatrix = cross->explicitmatrix;
300:   return(0);
301: }

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

306:    Not Collective

308:    Input Parameter:
309: .  svd  - singular value solver

311:    Output Parameter:
312: .  explicit - the mode flag

314:    Level: advanced

316: .seealso: SVDCrossSetExplicitMatrix()
317: @*/
318: PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
319: {

325:   PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
326:   return(0);
327: }

329: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
330: {
332:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

335:   PetscObjectReference((PetscObject)eps);
336:   EPSDestroy(&cross->eps);
337:   cross->eps = eps;
338:   cross->usereps = PETSC_TRUE;
339:   PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
340:   svd->state = SVD_STATE_INITIAL;
341:   return(0);
342: }

344: /*@
345:    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
346:    singular value solver.

348:    Collective on svd

350:    Input Parameters:
351: +  svd - singular value solver
352: -  eps - the eigensolver object

354:    Level: advanced

356: .seealso: SVDCrossGetEPS()
357: @*/
358: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
359: {

366:   PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
367:   return(0);
368: }

370: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
371: {
372:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

376:   if (!cross->eps) {
377:     EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
378:     PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
379:     EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
380:     EPSAppendOptionsPrefix(cross->eps,"svd_cross_");
381:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
382:     PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options);
383:     EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
384:     EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL);
385:   }
386:   *eps = cross->eps;
387:   return(0);
388: }

390: /*@
391:    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
392:    to the singular value solver.

394:    Not Collective

396:    Input Parameter:
397: .  svd - singular value solver

399:    Output Parameter:
400: .  eps - the eigensolver object

402:    Level: advanced

404: .seealso: SVDCrossSetEPS()
405: @*/
406: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
407: {

413:   PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
414:   return(0);
415: }

417: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
418: {
420:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
421:   PetscBool      isascii;

424:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
425:   if (isascii) {
426:     if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
427:     PetscViewerASCIIPrintf(viewer,"  %s matrix\n",cross->explicitmatrix?"explicit":"implicit");
428:     PetscViewerASCIIPushTab(viewer);
429:     EPSView(cross->eps,viewer);
430:     PetscViewerASCIIPopTab(viewer);
431:   }
432:   return(0);
433: }

435: PetscErrorCode SVDReset_Cross(SVD svd)
436: {
438:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

441:   EPSReset(cross->eps);
442:   MatDestroy(&cross->mat);
443:   VecDestroy(&cross->w);
444:   VecDestroy(&cross->diag);
445:   return(0);
446: }

448: PetscErrorCode SVDDestroy_Cross(SVD svd)
449: {
451:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

454:   EPSDestroy(&cross->eps);
455:   PetscFree(svd->data);
456:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
457:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
458:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL);
459:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL);
460:   return(0);
461: }

463: SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
464: {
466:   SVD_CROSS      *cross;

469:   PetscNewLog(svd,&cross);
470:   svd->data = (void*)cross;

472:   svd->ops->solve          = SVDSolve_Cross;
473:   svd->ops->setup          = SVDSetUp_Cross;
474:   svd->ops->setfromoptions = SVDSetFromOptions_Cross;
475:   svd->ops->destroy        = SVDDestroy_Cross;
476:   svd->ops->reset          = SVDReset_Cross;
477:   svd->ops->view           = SVDView_Cross;
478:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
479:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
480:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross);
481:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross);
482:   return(0);
483: }