Actual source code: cyclic.c

  1: /*                       

  3:    SLEPc singular value solver: "cyclic"

  5:    Method: Uses a Hermitian eigensolver for H(A) = [ 0  A ; A^T 0 ]

  7:    Last update: Jun 2007

  9:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 10:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 11:    Copyright (c) 2002-2009, 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 slepceps.h

 32: typedef struct {
 33:   PetscTruth explicitmatrix;
 34:   EPS        eps;
 35:   Mat        mat;
 36:   Vec        x1,x2,y1,y2;
 37: } SVD_CYCLIC;

 41: PetscErrorCode ShellMatMult_CYCLIC(Mat B,Vec x, Vec y)
 42: {
 44:   SVD            svd;
 45:   SVD_CYCLIC     *cyclic;
 46:   PetscScalar    *px,*py;
 47:   PetscInt       m;
 48: 
 50:   MatShellGetContext(B,(void**)&svd);
 51:   cyclic = (SVD_CYCLIC *)svd->data;
 52:   SVDMatGetLocalSize(svd,&m,PETSC_NULL);
 53:   VecGetArray(x,&px);
 54:   VecGetArray(y,&py);
 55:   VecPlaceArray(cyclic->x1,px);
 56:   VecPlaceArray(cyclic->x2,px+m);
 57:   VecPlaceArray(cyclic->y1,py);
 58:   VecPlaceArray(cyclic->y2,py+m);
 59:   SVDMatMult(svd,PETSC_FALSE,cyclic->x2,cyclic->y1);
 60:   SVDMatMult(svd,PETSC_TRUE,cyclic->x1,cyclic->y2);
 61:   VecResetArray(cyclic->x1);
 62:   VecResetArray(cyclic->x2);
 63:   VecResetArray(cyclic->y1);
 64:   VecResetArray(cyclic->y2);
 65:   VecRestoreArray(x,&px);
 66:   VecRestoreArray(y,&py);
 67:   return(0);
 68: }

 72: PetscErrorCode ShellMatGetDiagonal_CYCLIC(Mat B,Vec diag)
 73: {
 75: 
 77:   VecSet(diag,0.0);
 78:   return(0);
 79: }


 84: PetscErrorCode SVDSetUp_CYCLIC(SVD svd)
 85: {
 86:   PetscErrorCode    ierr;
 87:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC *)svd->data;
 88:   PetscInt          M,N,m,n,i,j,start,end,ncols,*pos,nloc;
 89:   const PetscInt    *cols;
 90:   const PetscScalar *vals;
 91:   PetscScalar       *pU;

 94: 
 95:   if (cyclic->mat) {
 96:     MatDestroy(cyclic->mat);
 97:   }
 98:   if (cyclic->x1) {
 99:     VecDestroy(cyclic->x1);
100:     VecDestroy(cyclic->x2);
101:     VecDestroy(cyclic->y1);
102:     VecDestroy(cyclic->y2);
103:   }

105:   SVDMatGetSize(svd,&M,&N);
106:   SVDMatGetLocalSize(svd,&m,&n);
107:   if (cyclic->explicitmatrix) {
108:     cyclic->x1 = cyclic->x2 = cyclic->y1 = cyclic->y2 = PETSC_NULL;
109:     MatCreate(((PetscObject)svd)->comm,&cyclic->mat);
110:     MatSetSizes(cyclic->mat,m+n,m+n,M+N,M+N);
111:     MatSetFromOptions(cyclic->mat);
112:     if (svd->AT) {
113:       MatGetOwnershipRange(svd->AT,&start,&end);
114:       for (i=start;i<end;i++) {
115:         MatGetRow(svd->AT,i,&ncols,&cols,&vals);
116:         j = i + M;
117:         MatSetValues(cyclic->mat,1,&j,ncols,cols,vals,INSERT_VALUES);
118:         MatSetValues(cyclic->mat,ncols,cols,1,&j,vals,INSERT_VALUES);
119:         MatRestoreRow(svd->AT,i,&ncols,&cols,&vals);
120:       }
121:     } else {
122:       PetscMalloc(sizeof(PetscInt)*n,&pos);
123:       MatGetOwnershipRange(svd->A,&start,&end);
124:       for (i=start;i<end;i++) {
125:         MatGetRow(svd->A,i,&ncols,&cols,&vals);
126:         for (j=0;j<ncols;j++)
127:           pos[j] = cols[j] + M;
128:         MatSetValues(cyclic->mat,1,&i,ncols,pos,vals,INSERT_VALUES);
129:         MatSetValues(cyclic->mat,ncols,pos,1,&i,vals,INSERT_VALUES);
130:         MatRestoreRow(svd->A,i,&ncols,&cols,&vals);
131:       }
132:       PetscFree(pos);
133:     }
134:     MatAssemblyBegin(cyclic->mat,MAT_FINAL_ASSEMBLY);
135:     MatAssemblyEnd(cyclic->mat,MAT_FINAL_ASSEMBLY);
136:   } else {
137:     VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->x1);
138:     VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->x2);
139:     VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->y1);
140:     VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->y2);
141:     MatCreateShell(((PetscObject)svd)->comm,m+n,m+n,M+N,M+N,svd,&cyclic->mat);
142:     MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))ShellMatMult_CYCLIC);
143:     MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_CYCLIC);
144:   }

146:   EPSSetOperators(cyclic->eps,cyclic->mat,PETSC_NULL);
147:   EPSSetProblemType(cyclic->eps,EPS_HEP);
148:   EPSSetWhichEigenpairs(cyclic->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_MAGNITUDE);
149:   EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
150:   EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
151:   EPSSetUp(cyclic->eps);
152:   EPSGetDimensions(cyclic->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
153:   EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);

155:   if (svd->ncv != svd->n) {
156:     if (svd->U) {
157:       VecGetArray(svd->U[0],&pU);
158:       for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]);  }
159:       PetscFree(pU);
160:       PetscFree(svd->U);
161:     }
162:     PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
163:     SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
164:     PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
165:     for (i=0;i<svd->ncv;i++) {
166:       VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
167:     }
168:   }

170:   return(0);
171: }

175: PetscErrorCode SVDSolve_CYCLIC(SVD svd)
176: {
178:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC *)svd->data;
179:   PetscInt       i,j,M,m,idx,start,end;
180:   PetscScalar    sigma,*px;
181:   Vec            x;
182:   IS             isU,isV;
183:   VecScatter     vsU,vsV;
184: 
186:   EPSSolve(cyclic->eps);
187:   EPSGetConverged(cyclic->eps,&svd->nconv);
188:   EPSGetIterationNumber(cyclic->eps,&svd->its);
189:   EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);

191:   MatGetVecs(cyclic->mat,&x,PETSC_NULL);
192:   if (cyclic->explicitmatrix) {
193:     EPSGetOperationCounters(cyclic->eps,&svd->matvecs,PETSC_NULL,PETSC_NULL);
194:     SVDMatGetSize(svd,&M,PETSC_NULL);
195:     VecGetOwnershipRange(svd->U[0],&start,&end);
196:     ISCreateBlock(((PetscObject)svd)->comm,end-start,1,&start,&isU);
197:     VecScatterCreate(x,isU,svd->U[0],PETSC_NULL,&vsU);

199:     VecGetOwnershipRange(svd->V[0],&start,&end);
200:     idx = start + M;
201:     ISCreateBlock(((PetscObject)svd)->comm,end-start,1,&idx,&isV);
202:     VecScatterCreate(x,isV,svd->V[0],PETSC_NULL,&vsV);

204:     for (i=0,j=0;i<svd->nconv;i++) {
205:       EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
206:       if (PetscRealPart(sigma) > 0.0) {
207:         svd->sigma[j] = PetscRealPart(sigma);
208:         VecScatterBegin(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
209:         VecScatterBegin(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
210:         VecScatterEnd(vsU,x,svd->U[j],INSERT_VALUES,SCATTER_FORWARD);
211:         VecScatterEnd(vsV,x,svd->V[j],INSERT_VALUES,SCATTER_FORWARD);
212:         VecScale(svd->U[j],1.0/sqrt(2.0));
213:         VecScale(svd->V[j],1.0/sqrt(2.0));
214:         j++;
215:       }
216:     }
217: 
218:     ISDestroy(isU);
219:     VecScatterDestroy(vsU);
220:     ISDestroy(isV);
221:     VecScatterDestroy(vsV);
222:   } else {
223:     SVDMatGetLocalSize(svd,&m,PETSC_NULL);
224:     for (i=0,j=0;i<svd->nconv;i++) {
225:       EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
226:       if (PetscRealPart(sigma) > 0.0) {
227:         svd->sigma[j] = PetscRealPart(sigma);
228:         VecGetArray(x,&px);
229:         VecPlaceArray(cyclic->x1,px);
230:         VecPlaceArray(cyclic->x2,px+m);
231: 
232:         VecCopy(cyclic->x1,svd->U[j]);
233:         VecScale(svd->U[j],1.0/sqrt(2.0));

235:         VecCopy(cyclic->x2,svd->V[j]);
236:         VecScale(svd->V[j],1.0/sqrt(2.0));
237: 
238:         VecResetArray(cyclic->x1);
239:         VecResetArray(cyclic->x2);
240:         VecRestoreArray(x,&px);
241:         j++;
242:       }
243:     }
244:   }
245:   svd->nconv = j;

247:   VecDestroy(x);
248:   return(0);
249: }

253: PetscErrorCode SVDMonitor_CYCLIC(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
254: {
255:   PetscInt   i,j;
256:   SVD        svd = (SVD)ctx;

259:   nconv = 0;
260:   for (i=0,j=0;i<nest;i++) {
261:     if (PetscRealPart(eigr[i]) > 0.0) {
262:       svd->sigma[j] = PetscRealPart(eigr[i]);
263:       svd->errest[j] = errest[i];
264:       if (errest[i] < svd->tol) nconv++;
265:       j++;
266:     }
267:   }
268:   nest = j;
269:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
270:   return(0);
271: }

275: PetscErrorCode SVDSetFromOptions_CYCLIC(SVD svd)
276: {
278:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC *)svd->data;
279:   ST             st;

282:   PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"CYCLIC Singular Value Solver Options","SVD");
283:   PetscOptionsTruth("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",PETSC_FALSE,&cyclic->explicitmatrix,PETSC_NULL);
284:   if (cyclic->explicitmatrix) {
285:     /* don't build the transpose */
286:     if (svd->transmode == PETSC_DECIDE)
287:       svd->transmode = SVD_TRANSPOSE_IMPLICIT;
288:   } else {
289:     /* use as default an ST with shell matrix and Jacobi */
290:     EPSGetST(cyclic->eps,&st);
291:     STSetMatMode(st,STMATMODE_SHELL);
292:   }
293:   PetscOptionsEnd();
294:   EPSSetFromOptions(cyclic->eps);
295:   PetscOptionsTail();
296:   return(0);
297: }

302: PetscErrorCode SVDCyclicSetExplicitMatrix_CYCLIC(SVD svd,PetscTruth explicitmatrix)
303: {
304:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

307:   cyclic->explicitmatrix = explicitmatrix;
308:   return(0);
309: }

314: /*@
315:    SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator 
316:    H(A) = [ 0  A ; A^T 0 ] must be computed explicitly.

318:    Collective on SVD

320:    Input Parameters:
321: +  svd      - singular value solver
322: -  explicit - boolean flag indicating if H(A) is built explicitly

324:    Options Database Key:
325: .  -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag

327:    Level: advanced

329: .seealso: SVDCyclicGetExplicitMatrix()
330: @*/
331: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscTruth explicitmatrix)
332: {
333:   PetscErrorCode ierr, (*f)(SVD,PetscTruth);

337:   PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",(void (**)())&f);
338:   if (f) {
339:     (*f)(svd,explicitmatrix);
340:   }
341:   return(0);
342: }

347: PetscErrorCode SVDCyclicGetExplicitMatrix_CYCLIC(SVD svd,PetscTruth *explicitmatrix)
348: {
349:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

353:   *explicitmatrix = cyclic->explicitmatrix;
354:   return(0);
355: }

360: /*@C
361:    SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly

363:    Not collective

365:    Input Parameter:
366: .  svd  - singular value solver

368:    Output Parameter:
369: .  explicit - the mode flag

371:    Level: advanced

373: .seealso: SVDCyclicSetExplicitMatrix()
374: @*/
375: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscTruth *explicitmatrix)
376: {
377:   PetscErrorCode ierr, (*f)(SVD,PetscTruth*);

381:   PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",(void (**)())&f);
382:   if (f) {
383:     (*f)(svd,explicitmatrix);
384:   }
385:   return(0);
386: }

391: PetscErrorCode SVDCyclicSetEPS_CYCLIC(SVD svd,EPS eps)
392: {
393:   PetscErrorCode  ierr;
394:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

399:   PetscObjectReference((PetscObject)eps);
400:   EPSDestroy(cyclic->eps);
401:   cyclic->eps = eps;
402:   svd->setupcalled = 0;
403:   return(0);
404: }

409: /*@
410:    SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
411:    singular value solver. 

413:    Collective on SVD

415:    Input Parameters:
416: +  svd - singular value solver
417: -  eps - the eigensolver object

419:    Level: advanced

421: .seealso: SVDCyclicGetEPS()
422: @*/
423: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
424: {
425:   PetscErrorCode ierr, (*f)(SVD,EPS eps);

429:   PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicSetEPS_C",(void (**)())&f);
430:   if (f) {
431:     (*f)(svd,eps);
432:   }
433:   return(0);
434: }

439: PetscErrorCode SVDCyclicGetEPS_CYCLIC(SVD svd,EPS *eps)
440: {
441:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

445:   *eps = cyclic->eps;
446:   return(0);
447: }

452: /*@C
453:    SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
454:    to the singular value solver.

456:    Not Collective

458:    Input Parameter:
459: .  svd - singular value solver

461:    Output Parameter:
462: .  eps - the eigensolver object

464:    Level: advanced

466: .seealso: SVDCyclicSetEPS()
467: @*/
468: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
469: {
470:   PetscErrorCode ierr, (*f)(SVD,EPS *eps);

474:   PetscObjectQueryFunction((PetscObject)svd,"SVDCyclicGetEPS_C",(void (**)())&f);
475:   if (f) {
476:     (*f)(svd,eps);
477:   }
478:   return(0);
479: }

483: PetscErrorCode SVDView_CYCLIC(SVD svd,PetscViewer viewer)
484: {
485:   PetscErrorCode  ierr;
486:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

489:   if (cyclic->explicitmatrix) {
490:     PetscViewerASCIIPrintf(viewer,"cyclic matrix: explicit\n");
491:   } else {
492:     PetscViewerASCIIPrintf(viewer,"cyclic matrix: implicit\n");
493:   }
494:   EPSView(cyclic->eps,viewer);
495:   return(0);
496: }

500: PetscErrorCode SVDDestroy_CYCLIC(SVD svd)
501: {
502:   PetscErrorCode  ierr;
503:   SVD_CYCLIC *cyclic = (SVD_CYCLIC *)svd->data;

506:   EPSDestroy(cyclic->eps);
507:   if (cyclic->mat) { MatDestroy(cyclic->mat); }
508:   if (cyclic->x1) {
509:     VecDestroy(cyclic->x1);
510:     VecDestroy(cyclic->x2);
511:     VecDestroy(cyclic->y1);
512:     VecDestroy(cyclic->y2);
513:   }
514:   PetscFree(svd->data);
515:   return(0);
516: }

521: PetscErrorCode SVDCreate_CYCLIC(SVD svd)
522: {
523:   PetscErrorCode  ierr;
524:   SVD_CYCLIC *cyclic;
525: 
527:   PetscNew(SVD_CYCLIC,&cyclic);
528:   PetscLogObjectMemory(svd,sizeof(SVD_CYCLIC));
529:   svd->data                      = (void *)cyclic;
530:   svd->ops->solve                = SVDSolve_CYCLIC;
531:   svd->ops->setup                = SVDSetUp_CYCLIC;
532:   svd->ops->setfromoptions       = SVDSetFromOptions_CYCLIC;
533:   svd->ops->destroy              = SVDDestroy_CYCLIC;
534:   svd->ops->view                 = SVDView_CYCLIC;
535:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","SVDCyclicSetEPS_CYCLIC",SVDCyclicSetEPS_CYCLIC);
536:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","SVDCyclicGetEPS_CYCLIC",SVDCyclicGetEPS_CYCLIC);
537:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","SVDCyclicSetExplicitMatrix_CYCLIC",SVDCyclicSetExplicitMatrix_CYCLIC);
538:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","SVDCyclicGetExplicitMatrix_CYCLIC",SVDCyclicGetExplicitMatrix_CYCLIC);

540:   EPSCreate(((PetscObject)svd)->comm,&cyclic->eps);
541:   EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
542:   EPSAppendOptionsPrefix(cyclic->eps,"svd_");
543:   PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
544:   PetscLogObjectParent(svd,cyclic->eps);
545:   EPSSetIP(cyclic->eps,svd->ip);
546:   EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
547:   EPSMonitorSet(cyclic->eps,SVDMonitor_CYCLIC,svd,PETSC_NULL);
548:   cyclic->explicitmatrix = PETSC_FALSE;
549:   cyclic->mat = PETSC_NULL;
550:   cyclic->x1 = PETSC_NULL;
551:   cyclic->x2 = PETSC_NULL;
552:   cyclic->y1 = PETSC_NULL;
553:   cyclic->y2 = PETSC_NULL;
554:   return(0);
555: }