Actual source code: cross.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*

  3:    SLEPc singular value solver: "cross"

  5:    Method: Uses a Hermitian eigensolver for A^T*A

  7:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  8:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  9:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 11:    This file is part of SLEPc.

 13:    SLEPc is free software: you can redistribute it and/or modify it under  the
 14:    terms of version 3 of the GNU Lesser General Public License as published by
 15:    the Free Software Foundation.

 17:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 18:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 19:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 20:    more details.

 22:    You  should have received a copy of the GNU Lesser General  Public  License
 23:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 24:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25: */

 27: #include <slepc/private/svdimpl.h>                /*I "slepcsvd.h" I*/
 28: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/

 30: typedef struct {
 31:   EPS       eps;
 32:   Mat       mat;
 33:   Vec       w,diag;
 34: } SVD_CROSS;

 38: static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
 39: {
 41:   SVD            svd;
 42:   SVD_CROSS      *cross;

 45:   MatShellGetContext(B,(void**)&svd);
 46:   cross = (SVD_CROSS*)svd->data;
 47:   SVDMatMult(svd,PETSC_FALSE,x,cross->w);
 48:   SVDMatMult(svd,PETSC_TRUE,cross->w,y);
 49:   return(0);
 50: }

 54: static PetscErrorCode MatCreateVecs_Cross(Mat B,Vec *right,Vec *left)
 55: {
 57:   SVD            svd;

 60:   MatShellGetContext(B,(void**)&svd);
 61:   if (right) {
 62:     SVDMatCreateVecs(svd,right,NULL);
 63:     if (left) { VecDuplicate(*right,left); }
 64:   } else {
 65:     SVDMatCreateVecs(svd,left,NULL);
 66:   }
 67:   return(0);
 68: }

 72: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
 73: {
 74:   PetscErrorCode    ierr;
 75:   SVD               svd;
 76:   SVD_CROSS         *cross;
 77:   PetscMPIInt       len;
 78:   PetscInt          N,n,i,j,start,end,ncols;
 79:   PetscScalar       *work1,*work2,*diag;
 80:   const PetscInt    *cols;
 81:   const PetscScalar *vals;

 84:   MatShellGetContext(B,(void**)&svd);
 85:   cross = (SVD_CROSS*)svd->data;
 86:   if (!cross->diag) {
 87:     /* compute diagonal from rows and store in cross->diag */
 88:     VecDuplicate(d,&cross->diag);
 89:     SVDMatGetSize(svd,NULL,&N);
 90:     SVDMatGetLocalSize(svd,NULL,&n);
 91:     PetscMalloc2(N,&work1,N,&work2);
 92:     for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
 93:     if (svd->AT) {
 94:       MatGetOwnershipRange(svd->AT,&start,&end);
 95:       for (i=start;i<end;i++) {
 96:         MatGetRow(svd->AT,i,&ncols,NULL,&vals);
 97:         for (j=0;j<ncols;j++)
 98:           work1[i] += vals[j]*vals[j];
 99:         MatRestoreRow(svd->AT,i,&ncols,NULL,&vals);
100:       }
101:     } else {
102:       MatGetOwnershipRange(svd->A,&start,&end);
103:       for (i=start;i<end;i++) {
104:         MatGetRow(svd->A,i,&ncols,&cols,&vals);
105:         for (j=0;j<ncols;j++)
106:           work1[cols[j]] += vals[j]*vals[j];
107:         MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
108:       }
109:     }
110:     PetscMPIIntCast(N,&len);
111:     MPI_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)svd));
112:     VecGetOwnershipRange(cross->diag,&start,&end);
113:     VecGetArray(cross->diag,&diag);
114:     for (i=start;i<end;i++) diag[i-start] = work2[i];
115:     VecRestoreArray(cross->diag,&diag);
116:     PetscFree2(work1,work2);
117:   }
118:   VecCopy(cross->diag,d);
119:   return(0);
120: }

124: PetscErrorCode SVDSetUp_Cross(SVD svd)
125: {
127:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
128:   PetscInt       n;
129:   PetscBool      trackall;

132:   if (!cross->mat) {
133:     SVDMatGetLocalSize(svd,NULL,&n);
134:     MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
135:     MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
136:     MatShellSetOperation(cross->mat,MATOP_GET_VECS,(void(*)(void))MatCreateVecs_Cross);
137:     MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
138:     SVDMatCreateVecs(svd,NULL,&cross->w);
139:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->mat);
140:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->w);
141:   }

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

177: PetscErrorCode SVDSolve_Cross(SVD svd)
178: {
180:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
181:   PetscInt       i;
182:   PetscScalar    sigma;
183:   Vec            v;

186:   EPSSolve(cross->eps);
187:   EPSGetConverged(cross->eps,&svd->nconv);
188:   EPSGetIterationNumber(cross->eps,&svd->its);
189:   EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
190:   for (i=0;i<svd->nconv;i++) {
191:     BVGetColumn(svd->V,i,&v);
192:     EPSGetEigenpair(cross->eps,i,&sigma,NULL,v,NULL);
193:     BVRestoreColumn(svd->V,i,&v);
194:     if (PetscRealPart(sigma)<0.0) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS");
195:     svd->sigma[i] = PetscSqrtReal(PetscRealPart(sigma));
196:   }
197:   return(0);
198: }

202: static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
203: {
204:   PetscInt       i;
205:   SVD            svd = (SVD)ctx;
206:   PetscScalar    er,ei;

210:   for (i=0;i<PetscMin(nest,svd->ncv);i++) {
211:     er = eigr[i]; ei = eigi[i];
212:     STBackTransform(eps->st,1,&er,&ei);
213:     svd->sigma[i] = PetscSqrtReal(PetscRealPart(er));
214:     svd->errest[i] = errest[i];
215:   }
216:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
217:   return(0);
218: }

222: PetscErrorCode SVDSetFromOptions_Cross(PetscOptionItems *PetscOptionsObject,SVD svd)
223: {
225:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

228:   PetscOptionsHead(PetscOptionsObject,"SVD Cross Options");
229:   if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
230:   EPSSetFromOptions(cross->eps);
231:   PetscOptionsTail();
232:   return(0);
233: }

237: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
238: {
240:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

243:   PetscObjectReference((PetscObject)eps);
244:   EPSDestroy(&cross->eps);
245:   cross->eps = eps;
246:   PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
247:   svd->state = SVD_STATE_INITIAL;
248:   return(0);
249: }

253: /*@
254:    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
255:    singular value solver.

257:    Collective on SVD

259:    Input Parameters:
260: +  svd - singular value solver
261: -  eps - the eigensolver object

263:    Level: advanced

265: .seealso: SVDCrossGetEPS()
266: @*/
267: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
268: {

275:   PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
276:   return(0);
277: }

281: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
282: {
283:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
284:   ST             st;

288:   if (!cross->eps) {
289:     EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
290:     EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
291:     EPSAppendOptionsPrefix(cross->eps,"svd_cross_");
292:     PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
293:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cross->eps);
294:     EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
295:     EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL);
296:     EPSGetST(cross->eps,&st);
297:     STSetMatMode(st,ST_MATMODE_SHELL);
298:   }
299:   *eps = cross->eps;
300:   return(0);
301: }

305: /*@
306:    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
307:    to the singular value solver.

309:    Not Collective

311:    Input Parameter:
312: .  svd - singular value solver

314:    Output Parameter:
315: .  eps - the eigensolver object

317:    Level: advanced

319: .seealso: SVDCrossSetEPS()
320: @*/
321: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
322: {

328:   PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
329:   return(0);
330: }

334: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
335: {
337:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
338:   PetscBool      isascii;

341:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
342:   if (isascii) {
343:     if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
344:     PetscViewerASCIIPushTab(viewer);
345:     EPSView(cross->eps,viewer);
346:     PetscViewerASCIIPopTab(viewer);
347:   }
348:   return(0);
349: }

353: PetscErrorCode SVDReset_Cross(SVD svd)
354: {
356:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

359:   if (cross->eps) { EPSReset(cross->eps); }
360:   MatDestroy(&cross->mat);
361:   VecDestroy(&cross->w);
362:   VecDestroy(&cross->diag);
363:   return(0);
364: }

368: PetscErrorCode SVDDestroy_Cross(SVD svd)
369: {
371:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

374:   EPSDestroy(&cross->eps);
375:   PetscFree(svd->data);
376:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
377:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
378:   return(0);
379: }

383: PETSC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
384: {
386:   SVD_CROSS      *cross;

389:   PetscNewLog(svd,&cross);
390:   svd->data = (void*)cross;

392:   svd->ops->solve          = SVDSolve_Cross;
393:   svd->ops->setup          = SVDSetUp_Cross;
394:   svd->ops->setfromoptions = SVDSetFromOptions_Cross;
395:   svd->ops->destroy        = SVDDestroy_Cross;
396:   svd->ops->reset          = SVDReset_Cross;
397:   svd->ops->view           = SVDView_Cross;
398:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
399:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
400:   return(0);
401: }