Actual source code: cross.c

  1: /*

  3:    SLEPc singular value solver: "cross"

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

  7:    Last update: Jun 2007

  9:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 10:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 11:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

 13:    This file is part of SLEPc.

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

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

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

 29: #include <slepc-private/svdimpl.h>                /*I "slepcsvd.h" I*/
 30: #include <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/

 32: typedef struct {
 33:   EPS       eps;
 34:   PetscBool setfromoptionscalled;
 35:   Mat       mat;
 36:   Vec       w,diag;
 37: } SVD_CROSS;

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

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

 57: static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
 58: {
 59:   PetscErrorCode    ierr;
 60:   SVD               svd;
 61:   SVD_CROSS         *cross;
 62:   PetscInt          N,n,i,j,start,end,ncols;
 63:   PetscScalar       *work1,*work2,*diag;
 64:   const PetscInt    *cols;
 65:   const PetscScalar *vals;

 68:   MatShellGetContext(B,(void**)&svd);
 69:   cross = (SVD_CROSS*)svd->data;
 70:   if (!cross->diag) {
 71:     /* compute diagonal from rows and store in cross->diag */
 72:     VecDuplicate(d,&cross->diag);
 73:     SVDMatGetSize(svd,NULL,&N);
 74:     SVDMatGetLocalSize(svd,NULL,&n);
 75:     PetscMalloc(sizeof(PetscScalar)*N,&work1);
 76:     PetscMalloc(sizeof(PetscScalar)*N,&work2);
 77:     for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
 78:     if (svd->AT) {
 79:       MatGetOwnershipRange(svd->AT,&start,&end);
 80:       for (i=start;i<end;i++) {
 81:         MatGetRow(svd->AT,i,&ncols,NULL,&vals);
 82:         for (j=0;j<ncols;j++)
 83:           work1[i] += vals[j]*vals[j];
 84:         MatRestoreRow(svd->AT,i,&ncols,NULL,&vals);
 85:       }
 86:     } else {
 87:       MatGetOwnershipRange(svd->A,&start,&end);
 88:       for (i=start;i<end;i++) {
 89:         MatGetRow(svd->A,i,&ncols,&cols,&vals);
 90:         for (j=0;j<ncols;j++)
 91:           work1[cols[j]] += vals[j]*vals[j];
 92:         MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
 93:       }
 94:     }
 95:     MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPI_SUM,PetscObjectComm((PetscObject)svd));
 96:     VecGetOwnershipRange(cross->diag,&start,&end);
 97:     VecGetArray(cross->diag,&diag);
 98:     for (i=start;i<end;i++)
 99:       diag[i-start] = work2[i];
100:     VecRestoreArray(cross->diag,&diag);
101:     PetscFree(work1);
102:     PetscFree(work2);
103:   }
104:   VecCopy(cross->diag,d);
105:   return(0);
106: }

110: PetscErrorCode SVDSetUp_Cross(SVD svd)
111: {
113:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
114:   PetscInt       n,i;
115:   PetscBool      trackall;

118:   if (!cross->mat) {
119:     SVDMatGetLocalSize(svd,NULL,&n);
120:     MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
121:     MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))MatMult_Cross);
122:     MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross);
123:     SVDMatGetVecs(svd,NULL,&cross->w);
124:     PetscLogObjectParent(svd,cross->mat);
125:     PetscLogObjectParent(svd,cross->w);
126:   }

128:   if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
129:   EPSSetOperators(cross->eps,cross->mat,NULL);
130:   EPSSetProblemType(cross->eps,EPS_HEP);
131:   EPSSetWhichEigenpairs(cross->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_REAL);
132:   EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd);
133:   EPSSetTolerances(cross->eps,svd->tol,svd->max_it);
134:   /* Transfer the trackall option from svd to eps */
135:   SVDGetTrackAll(svd,&trackall);
136:   EPSSetTrackAll(cross->eps,trackall);
137:   if (cross->setfromoptionscalled) {
138:     EPSSetFromOptions(cross->eps);
139:     cross->setfromoptionscalled = PETSC_FALSE;
140:   }
141:   EPSSetUp(cross->eps);
142:   EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd);
143:   EPSGetTolerances(cross->eps,&svd->tol,&svd->max_it);
144:   /* Transfer the initial space from svd to eps */
145:   if (svd->nini < 0) {
146:     EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS);
147:     for (i=0;i<-svd->nini;i++) {
148:       VecDestroy(&svd->IS[i]);
149:     }
150:     PetscFree(svd->IS);
151:     svd->nini = 0;
152:   }
153:   return(0);
154: }

158: PetscErrorCode SVDSolve_Cross(SVD svd)
159: {
161:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
162:   PetscInt       i;
163:   PetscScalar    sigma;

166:   EPSSolve(cross->eps);
167:   EPSGetConverged(cross->eps,&svd->nconv);
168:   EPSGetIterationNumber(cross->eps,&svd->its);
169:   EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
170:   for (i=0;i<svd->nconv;i++) {
171:     EPSGetEigenpair(cross->eps,i,&sigma,NULL,svd->V[i],NULL);
172:     if (PetscRealPart(sigma)<0.0) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Negative eigenvalue computed by EPS");
173:     svd->sigma[i] = PetscSqrtReal(PetscRealPart(sigma));
174:   }
175:   return(0);
176: }

180: static PetscErrorCode SVDMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
181: {
182:   PetscInt       i;
183:   SVD            svd = (SVD)ctx;
184:   PetscScalar    er,ei;

188:   for (i=0;i<PetscMin(nest,svd->ncv);i++) {
189:     er = eigr[i]; ei = eigi[i];
190:     STBackTransform(eps->st,1,&er,&ei);
191:     svd->sigma[i] = PetscSqrtReal(PetscRealPart(er));
192:     svd->errest[i] = errest[i];
193:   }
194:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
195:   return(0);
196: }

200: PetscErrorCode SVDSetFromOptions_Cross(SVD svd)
201: {
202:   SVD_CROSS *cross = (SVD_CROSS*)svd->data;

205:   cross->setfromoptionscalled = PETSC_TRUE;
206:   return(0);
207: }

211: static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
212: {
214:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

217:   PetscObjectReference((PetscObject)eps);
218:   EPSDestroy(&cross->eps);
219:   cross->eps = eps;
220:   PetscLogObjectParent(svd,cross->eps);
221:   svd->setupcalled = 0;
222:   return(0);
223: }

227: /*@
228:    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
229:    singular value solver.

231:    Collective on SVD

233:    Input Parameters:
234: +  svd - singular value solver
235: -  eps - the eigensolver object

237:    Level: advanced

239: .seealso: SVDCrossGetEPS()
240: @*/
241: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
242: {

249:   PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
250:   return(0);
251: }

255: static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
256: {
257:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
258:   ST             st;

262:   if (!cross->eps) {
263:     EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps);
264:     EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
265:     EPSAppendOptionsPrefix(cross->eps,"svd_");
266:     PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
267:     PetscLogObjectParent(svd,cross->eps);
268:     if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
269:     EPSSetIP(cross->eps,svd->ip);
270:     EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
271:     EPSMonitorSet(cross->eps,SVDMonitor_Cross,svd,NULL);
272:     EPSGetST(cross->eps,&st);
273:     STSetMatMode(st,ST_MATMODE_SHELL);
274:   }
275:   *eps = cross->eps;
276:   return(0);
277: }

281: /*@
282:    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
283:    to the singular value solver.

285:    Not Collective

287:    Input Parameter:
288: .  svd - singular value solver

290:    Output Parameter:
291: .  eps - the eigensolver object

293:    Level: advanced

295: .seealso: SVDCrossSetEPS()
296: @*/
297: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
298: {

304:   PetscTryMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
305:   return(0);
306: }

310: PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
311: {
313:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

316:   if (!cross->eps) { SVDCrossGetEPS(svd,&cross->eps); }
317:   PetscViewerASCIIPushTab(viewer);
318:   EPSView(cross->eps,viewer);
319:   PetscViewerASCIIPopTab(viewer);
320:   return(0);
321: }

325: PetscErrorCode SVDReset_Cross(SVD svd)
326: {
328:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

331:   if (cross->eps) { EPSReset(cross->eps); }
332:   MatDestroy(&cross->mat);
333:   VecDestroy(&cross->w);
334:   VecDestroy(&cross->diag);
335:   return(0);
336: }

340: PetscErrorCode SVDDestroy_Cross(SVD svd)
341: {
343:   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;

346:   EPSDestroy(&cross->eps);
347:   PetscFree(svd->data);
348:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL);
349:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL);
350:   return(0);
351: }

355: PETSC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
356: {
358:   SVD_CROSS      *cross;

361:   PetscNewLog(svd,SVD_CROSS,&cross);
362:   svd->data                = (void*)cross;
363:   svd->ops->solve          = SVDSolve_Cross;
364:   svd->ops->setup          = SVDSetUp_Cross;
365:   svd->ops->setfromoptions = SVDSetFromOptions_Cross;
366:   svd->ops->destroy        = SVDDestroy_Cross;
367:   svd->ops->reset          = SVDReset_Cross;
368:   svd->ops->view           = SVDView_Cross;
369:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross);
370:   PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross);
371:   return(0);
372: }