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-2010, Universidad Politecnica de Valencia, Spain

 13:    This file is part of SLEPc.
 14:       
 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 private/svdimpl.h
 30:  #include private/epsimpl.h
 31:  #include slepceps.h

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

 42: PetscErrorCode ShellMatMult_CROSS(Mat B,Vec x, Vec y)
 43: {
 45:   SVD            svd;
 46:   SVD_CROSS      *cross;
 47: 
 49:   MatShellGetContext(B,(void**)&svd);
 50:   cross = (SVD_CROSS *)svd->data;
 51:   SVDMatMult(svd,PETSC_FALSE,x,cross->w);
 52:   SVDMatMult(svd,PETSC_TRUE,cross->w,y);
 53:   return(0);
 54: }

 58: PetscErrorCode ShellMatGetDiagonal_CROSS(Mat B,Vec d)
 59: {
 60:   PetscErrorCode    ierr;
 61:   SVD               svd;
 62:   SVD_CROSS         *cross;
 63:   PetscInt          N,n,i,j,start,end,ncols;
 64:   PetscScalar       *work1,*work2,*diag;
 65:   const PetscInt    *cols;
 66:   const PetscScalar *vals;
 67: 
 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,PETSC_NULL,&N);
 75:     SVDMatGetLocalSize(svd,PETSC_NULL,&n);
 76:     PetscMalloc(sizeof(PetscScalar)*N,&work1);
 77:     PetscMalloc(sizeof(PetscScalar)*N,&work2);
 78:     for (i=0;i<n;i++) work1[i] = work2[i] = 0.0;
 79:     if (svd->AT) {
 80:       MatGetOwnershipRange(svd->AT,&start,&end);
 81:       for (i=start;i<end;i++) {
 82:         MatGetRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
 83:         for (j=0;j<ncols;j++)
 84:           work1[i] += vals[j]*vals[j];
 85:         MatRestoreRow(svd->AT,i,&ncols,PETSC_NULL,&vals);
 86:       }
 87:     } else {
 88:       MatGetOwnershipRange(svd->A,&start,&end);
 89:       for (i=start;i<end;i++) {
 90:         MatGetRow(svd->A,i,&ncols,&cols,&vals);
 91:         for (j=0;j<ncols;j++)
 92:           work1[cols[j]] += vals[j]*vals[j];
 93:         MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
 94:       }
 95:     }
 96:     MPI_Allreduce(work1,work2,N,MPIU_SCALAR,MPI_SUM,((PetscObject)svd)->comm);
 97:     VecGetOwnershipRange(cross->diag,&start,&end);
 98:     VecGetArray(cross->diag,&diag);
 99:     for (i=start;i<end;i++)
100:       diag[i-start] = work2[i];
101:     VecRestoreArray(cross->diag,&diag);
102:     PetscFree(work1);
103:     PetscFree(work2);
104:   }
105:   VecCopy(cross->diag,d);
106:   return(0);
107: }

111: PetscErrorCode SVDSetUp_CROSS(SVD svd)
112: {
113:   PetscErrorCode    ierr;
114:   SVD_CROSS         *cross = (SVD_CROSS *)svd->data;
115:   PetscInt          n,i;
116:   PetscTruth        trackall;

119:   if (cross->mat) {
120:      MatDestroy(cross->mat);
121:      VecDestroy(cross->w);
122:   }
123:   if (cross->diag) {
124:      VecDestroy(cross->diag);
125:   }
126: 
127:   SVDMatGetLocalSize(svd,PETSC_NULL,&n);
128:   MatCreateShell(((PetscObject)svd)->comm,n,n,PETSC_DETERMINE,PETSC_DETERMINE,svd,&cross->mat);
129:   MatShellSetOperation(cross->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CROSS);
130:   MatShellSetOperation(cross->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CROSS);
131:   SVDMatGetVecs(svd,PETSC_NULL,&cross->w);

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

162: PetscErrorCode SVDSolve_CROSS(SVD svd)
163: {
165:   SVD_CROSS      *cross = (SVD_CROSS *)svd->data;
166:   PetscInt       i;
167:   PetscScalar    sigma;
168: 
170:   EPSSolve(cross->eps);
171:   EPSGetConverged(cross->eps,&svd->nconv);
172:   EPSGetIterationNumber(cross->eps,&svd->its);
173:   EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason);
174:   for (i=0;i<svd->nconv;i++) {
175:     EPSGetEigenpair(cross->eps,i,&sigma,PETSC_NULL,svd->V[i],PETSC_NULL);
176:     svd->sigma[i] = sqrt(PetscRealPart(sigma));
177:   }
178:   return(0);
179: }

183: PetscErrorCode SVDMonitor_CROSS(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
184: {
185:   PetscInt       i;
186:   SVD            svd = (SVD)ctx;
187:   PetscScalar    er,ei;

191:   for (i=0;i<nest;i++) {
192:     er = eigr[i]; ei = eigi[i];
193:     STBackTransform(eps->OP, 1, &er, &ei);
194:     svd->sigma[i] = sqrt(PetscRealPart(er));
195:     svd->errest[i] = errest[i];
196:   }
197:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
198:   return(0);
199: }

203: PetscErrorCode SVDSetFromOptions_CROSS(SVD svd)
204: {
205:   SVD_CROSS      *cross = (SVD_CROSS *)svd->data;

208:   cross->setfromoptionscalled = PETSC_TRUE;
209:   return(0);
210: }

215: PetscErrorCode SVDCrossSetEPS_CROSS(SVD svd,EPS eps)
216: {
218:   SVD_CROSS      *cross = (SVD_CROSS *)svd->data;

223:   PetscObjectReference((PetscObject)eps);
224:   EPSDestroy(cross->eps);
225:   cross->eps = eps;
226:   svd->setupcalled = 0;
227:   return(0);
228: }

233: /*@
234:    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
235:    singular value solver. 

237:    Collective on SVD

239:    Input Parameters:
240: +  svd - singular value solver
241: -  eps - the eigensolver object

243:    Level: advanced

245: .seealso: SVDCrossGetEPS()
246: @*/
247: PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
248: {
249:   PetscErrorCode ierr, (*f)(SVD,EPS eps);

253:   PetscObjectQueryFunction((PetscObject)svd,"SVDCrossSetEPS_C",(void (**)())&f);
254:   if (f) {
255:     (*f)(svd,eps);
256:   }
257:   return(0);
258: }

263: PetscErrorCode SVDCrossGetEPS_CROSS(SVD svd,EPS *eps)
264: {
265:   SVD_CROSS *cross = (SVD_CROSS *)svd->data;

269:   *eps = cross->eps;
270:   return(0);
271: }

276: /*@
277:    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
278:    to the singular value solver.

280:    Not Collective

282:    Input Parameter:
283: .  svd - singular value solver

285:    Output Parameter:
286: .  eps - the eigensolver object

288:    Level: advanced

290: .seealso: SVDCrossSetEPS()
291: @*/
292: PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
293: {
294:   PetscErrorCode ierr, (*f)(SVD,EPS *eps);

298:   PetscObjectQueryFunction((PetscObject)svd,"SVDCrossGetEPS_C",(void (**)())&f);
299:   if (f) {
300:     (*f)(svd,eps);
301:   }
302:   return(0);
303: }

307: PetscErrorCode SVDView_CROSS(SVD svd,PetscViewer viewer)
308: {
310:   SVD_CROSS      *cross = (SVD_CROSS *)svd->data;

313:   EPSView(cross->eps,viewer);
314:   return(0);
315: }

319: PetscErrorCode SVDDestroy_CROSS(SVD svd)
320: {
322:   SVD_CROSS      *cross = (SVD_CROSS *)svd->data;

325:   EPSDestroy(cross->eps);
326:   if (cross->mat) {
327:     MatDestroy(cross->mat);
328:     VecDestroy(cross->w);
329:   }
330:   if (cross->diag) {
331:     VecDestroy(cross->diag);
332:   }
333:   PetscFree(svd->data);
334:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","",PETSC_NULL);
335:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","",PETSC_NULL);
336:   return(0);
337: }

342: PetscErrorCode SVDCreate_CROSS(SVD svd)
343: {
345:   SVD_CROSS      *cross;
346:   ST             st;
347: 
349:   PetscNew(SVD_CROSS,&cross);
350:   PetscLogObjectMemory(svd,sizeof(SVD_CROSS));
351:   svd->data                = (void *)cross;
352:   svd->ops->solve          = SVDSolve_CROSS;
353:   svd->ops->setup          = SVDSetUp_CROSS;
354:   svd->ops->setfromoptions = SVDSetFromOptions_CROSS;
355:   svd->ops->destroy        = SVDDestroy_CROSS;
356:   svd->ops->view           = SVDView_CROSS;
357:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossSetEPS_C","SVDCrossSetEPS_CROSS",SVDCrossSetEPS_CROSS);
358:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCrossGetEPS_C","SVDCrossGetEPS_CROSS",SVDCrossGetEPS_CROSS);

360:   EPSCreate(((PetscObject)svd)->comm,&cross->eps);
361:   EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix);
362:   EPSAppendOptionsPrefix(cross->eps,"svd_");
363:   PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1);
364:   PetscLogObjectParent(svd,cross->eps);
365:   EPSSetIP(cross->eps,svd->ip);
366:   EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL);
367:   EPSMonitorSet(cross->eps,SVDMonitor_CROSS,svd,PETSC_NULL);
368:   EPSGetST(cross->eps,&st);
369:   STSetMatMode(st,ST_MATMODE_SHELL);
370:   cross->mat = PETSC_NULL;
371:   cross->w = PETSC_NULL;
372:   cross->diag = PETSC_NULL;
373:   cross->setfromoptionscalled = PETSC_FALSE;
374:   return(0);
375: }