Actual source code: cyclic.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*

  3:    SLEPc singular value solver: "cyclic"

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

  7:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  8:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  9:    Copyright (c) 2002-2014, 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:   PetscBool explicitmatrix;
 32:   EPS       eps;
 33:   Mat       mat;
 34:   Vec       x1,x2,y1,y2;
 35: } SVD_CYCLIC;

 39: static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
 40: {
 41:   PetscErrorCode    ierr;
 42:   SVD               svd;
 43:   SVD_CYCLIC        *cyclic;
 44:   const PetscScalar *px;
 45:   PetscScalar       *py;
 46:   PetscInt          m;

 49:   MatShellGetContext(B,(void**)&svd);
 50:   cyclic = (SVD_CYCLIC*)svd->data;
 51:   SVDMatGetLocalSize(svd,&m,NULL);
 52:   VecGetArrayRead(x,&px);
 53:   VecGetArray(y,&py);
 54:   VecPlaceArray(cyclic->x1,px);
 55:   VecPlaceArray(cyclic->x2,px+m);
 56:   VecPlaceArray(cyclic->y1,py);
 57:   VecPlaceArray(cyclic->y2,py+m);
 58:   SVDMatMult(svd,PETSC_FALSE,cyclic->x2,cyclic->y1);
 59:   SVDMatMult(svd,PETSC_TRUE,cyclic->x1,cyclic->y2);
 60:   VecResetArray(cyclic->x1);
 61:   VecResetArray(cyclic->x2);
 62:   VecResetArray(cyclic->y1);
 63:   VecResetArray(cyclic->y2);
 64:   VecRestoreArrayRead(x,&px);
 65:   VecRestoreArray(y,&py);
 66:   return(0);
 67: }

 71: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
 72: {

 76:   VecSet(diag,0.0);
 77:   return(0);
 78: }

 82: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
 83: {
 84:   PetscErrorCode    ierr;
 85:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC*)svd->data;
 86:   PetscInt          M,N,m,n,i,isl,Istart,Iend;
 87:   const PetscScalar *isa;
 88:   PetscScalar       *va;
 89:   PetscBool         trackall;
 90:   Vec               v;
 91:   Mat               Zm,Zn;

 94:   SVDMatGetSize(svd,&M,&N);
 95:   SVDMatGetLocalSize(svd,&m,&n);
 96:   if (!cyclic->mat) {
 97:     if (cyclic->explicitmatrix) {
 98:       if (!svd->AT) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
 99:       MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
100:       MatSetSizes(Zm,m,m,M,M);
101:       MatSetFromOptions(Zm);
102:       MatSetUp(Zm);
103:       MatGetOwnershipRange(Zm,&Istart,&Iend);
104:       for (i=Istart;i<Iend;i++) {
105:         MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
106:       }
107:       MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
108:       MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
109:       MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
110:       MatSetSizes(Zn,n,n,N,N);
111:       MatSetFromOptions(Zn);
112:       MatSetUp(Zn);
113:       MatGetOwnershipRange(Zn,&Istart,&Iend);
114:       for (i=Istart;i<Iend;i++) {
115:         MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
116:       }
117:       MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
118:       MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
119:       SlepcMatTile(1.0,Zm,1.0,svd->A,1.0,svd->AT,1.0,Zn,&cyclic->mat);
120:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->mat);
121:       MatDestroy(&Zm);
122:       MatDestroy(&Zn);
123:     } else {
124:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&cyclic->x1);
125:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&cyclic->x2);
126:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&cyclic->y1);
127:       VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&cyclic->y2);
128:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->x1);
129:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->x2);
130:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->y1);
131:       PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->y2);
132:       MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,svd,&cyclic->mat);
133:       MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))MatMult_Cyclic);
134:       MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic);
135:     }
136:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->mat);
137:   }

139:   if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
140:   EPSSetOperators(cyclic->eps,cyclic->mat,NULL);
141:   EPSSetProblemType(cyclic->eps,EPS_HEP);
142:   if (svd->which == SVD_LARGEST) {
143:     EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
144:   } else {
145:     EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
146:     EPSSetTarget(cyclic->eps,0.0);
147:   }
148:   EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv?svd->ncv:PETSC_DEFAULT,svd->mpd?svd->mpd:PETSC_DEFAULT);
149:   EPSSetTolerances(cyclic->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it?svd->max_it:PETSC_DEFAULT);
150:   /* Transfer the trackall option from svd to eps */
151:   SVDGetTrackAll(svd,&trackall);
152:   EPSSetTrackAll(cyclic->eps,trackall);
153:   /* Transfer the initial subspace from svd to eps */
154:   if (svd->nini<0 || svd->ninil<0) {
155:     for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
156:       MatGetVecs(cyclic->mat,&v,NULL);
157:       VecGetArray(v,&va);
158:       if (i<-svd->ninil) {
159:         VecGetSize(svd->ISL[i],&isl);
160:         if (isl!=m) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
161:         VecGetArrayRead(svd->ISL[i],&isa);
162:         PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
163:         VecRestoreArrayRead(svd->IS[i],&isa);
164:       } else {
165:         PetscMemzero(&va,sizeof(PetscScalar)*m);
166:       }
167:       if (i<-svd->nini) {
168:         VecGetSize(svd->IS[i],&isl);
169:         if (isl!=n) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
170:         VecGetArrayRead(svd->IS[i],&isa);
171:         PetscMemcpy(va+m,isa,sizeof(PetscScalar)*n);
172:         VecRestoreArrayRead(svd->IS[i],&isa);
173:       } else {
174:         PetscMemzero(va+m,sizeof(PetscScalar)*n);
175:       }
176:       VecRestoreArray(v,&va);
177:       VecDestroy(&svd->IS[i]);
178:       svd->IS[i] = v;
179:     }
180:     svd->nini = PetscMin(svd->nini,svd->ninil);
181:     EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
182:     SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
183:     SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
184:   }
185:   EPSSetUp(cyclic->eps);
186:   EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
187:   svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
188:   EPSGetTolerances(cyclic->eps,NULL,&svd->max_it);
189:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

191:   svd->leftbasis = PETSC_TRUE;
192:   SVDAllocateSolution(svd,0);
193:   return(0);
194: }

198: PetscErrorCode SVDSolve_Cyclic(SVD svd)
199: {
200:   PetscErrorCode    ierr;
201:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC*)svd->data;
202:   PetscInt          i,j,M,N,m,n;
203:   PetscScalar       sigma;
204:   const PetscScalar *px;
205:   Vec               x,x1,x2;

208:   EPSSolve(cyclic->eps);
209:   EPSGetConverged(cyclic->eps,&svd->nconv);
210:   EPSGetIterationNumber(cyclic->eps,&svd->its);
211:   EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);

213:   MatGetVecs(cyclic->mat,&x,NULL);
214:   SVDMatGetSize(svd,&M,&N);
215:   SVDMatGetLocalSize(svd,&m,&n);
216:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,m,M,NULL,&x1);
217:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)svd),1,n,N,NULL,&x2);
218:   for (i=0,j=0;i<svd->nconv;i++) {
219:     EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
220:     if (PetscRealPart(sigma) > 0.0) {
221:       svd->sigma[j] = PetscRealPart(sigma);
222:       VecGetArrayRead(x,&px);
223:       VecPlaceArray(x1,px);
224:       VecPlaceArray(x2,px+m);
225:       BVInsertVec(svd->U,j,x1);
226:       BVScaleColumn(svd->U,j,1.0/PetscSqrtReal(2.0));
227:       BVInsertVec(svd->V,j,x2);
228:       BVScaleColumn(svd->V,j,1.0/PetscSqrtReal(2.0));
229:       VecResetArray(x1);
230:       VecResetArray(x2);
231:       VecRestoreArrayRead(x,&px);
232:       j++;
233:     }
234:   }
235:   svd->nconv = j;

237:   VecDestroy(&x);
238:   VecDestroy(&x1);
239:   VecDestroy(&x2);
240:   return(0);
241: }

245: static PetscErrorCode SVDMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
246: {
247:   PetscInt       i,j;
248:   SVD            svd = (SVD)ctx;
249:   PetscScalar    er,ei;

253:   nconv = 0;
254:   for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
255:     er = eigr[i]; ei = eigi[i];
256:     STBackTransform(eps->st,1,&er,&ei);
257:     if (PetscRealPart(er) > 0.0) {
258:       svd->sigma[j] = PetscRealPart(er);
259:       svd->errest[j] = errest[i];
260:       if (errest[i] < svd->tol) nconv++;
261:       j++;
262:     }
263:   }
264:   nest = j;
265:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
266:   return(0);
267: }

271: PetscErrorCode SVDSetFromOptions_Cyclic(SVD svd)
272: {
274:   PetscBool      set,val;
275:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;
276:   ST             st;

279:   if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
280:   EPSSetFromOptions(cyclic->eps);
281:   PetscOptionsHead("SVD Cyclic Options");
282:   PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
283:   if (set) {
284:     SVDCyclicSetExplicitMatrix(svd,val);
285:   }
286:   if (!cyclic->explicitmatrix) {
287:     /* use as default an ST with shell matrix and Jacobi */
288:     if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
289:     EPSGetST(cyclic->eps,&st);
290:     STSetMatMode(st,ST_MATMODE_SHELL);
291:   }
292:   PetscOptionsTail();
293:   return(0);
294: }

298: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
299: {
300:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

303:   cyclic->explicitmatrix = explicitmatrix;
304:   return(0);
305: }

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

313:    Logically Collective on SVD

315:    Input Parameters:
316: +  svd      - singular value solver
317: -  explicit - boolean flag indicating if H(A) is built explicitly

319:    Options Database Key:
320: .  -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag

322:    Level: advanced

324: .seealso: SVDCyclicGetExplicitMatrix()
325: @*/
326: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
327: {

333:   PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
334:   return(0);
335: }

339: static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
340: {
341:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

344:   *explicitmatrix = cyclic->explicitmatrix;
345:   return(0);
346: }

350: /*@
351:    SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly

353:    Not Collective

355:    Input Parameter:
356: .  svd  - singular value solver

358:    Output Parameter:
359: .  explicit - the mode flag

361:    Level: advanced

363: .seealso: SVDCyclicSetExplicitMatrix()
364: @*/
365: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
366: {

372:   PetscTryMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
373:   return(0);
374: }

378: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
379: {
380:   PetscErrorCode  ierr;
381:   SVD_CYCLIC      *cyclic = (SVD_CYCLIC*)svd->data;

384:   PetscObjectReference((PetscObject)eps);
385:   EPSDestroy(&cyclic->eps);
386:   cyclic->eps = eps;
387:   PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
388:   svd->setupcalled = 0;
389:   return(0);
390: }

394: /*@
395:    SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
396:    singular value solver.

398:    Collective on SVD

400:    Input Parameters:
401: +  svd - singular value solver
402: -  eps - the eigensolver object

404:    Level: advanced

406: .seealso: SVDCyclicGetEPS()
407: @*/
408: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
409: {

416:   PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
417:   return(0);
418: }

422: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
423: {
425:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

428:   if (!cyclic->eps) {
429:     EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
430:     EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
431:     EPSAppendOptionsPrefix(cyclic->eps,"svd_");
432:     PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
433:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
434:     EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
435:     EPSMonitorSet(cyclic->eps,SVDMonitor_Cyclic,svd,NULL);
436:   }
437:   *eps = cyclic->eps;
438:   return(0);
439: }

443: /*@
444:    SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
445:    to the singular value solver.

447:    Not Collective

449:    Input Parameter:
450: .  svd - singular value solver

452:    Output Parameter:
453: .  eps - the eigensolver object

455:    Level: advanced

457: .seealso: SVDCyclicSetEPS()
458: @*/
459: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
460: {

466:   PetscTryMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
467:   return(0);
468: }

472: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
473: {
475:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

478:   if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
479:   PetscViewerASCIIPrintf(viewer,"  Cyclic: %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
480:   PetscViewerASCIIPushTab(viewer);
481:   EPSView(cyclic->eps,viewer);
482:   PetscViewerASCIIPopTab(viewer);
483:   return(0);
484: }

488: PetscErrorCode SVDReset_Cyclic(SVD svd)
489: {
491:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

494:   if (!cyclic->eps) { EPSReset(cyclic->eps); }
495:   MatDestroy(&cyclic->mat);
496:   VecDestroy(&cyclic->x1);
497:   VecDestroy(&cyclic->x2);
498:   VecDestroy(&cyclic->y1);
499:   VecDestroy(&cyclic->y2);
500:   return(0);
501: }

505: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
506: {
508:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

511:   EPSDestroy(&cyclic->eps);
512:   PetscFree(svd->data);
513:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
514:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
515:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
516:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
517:   return(0);
518: }

522: PETSC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
523: {
525:   SVD_CYCLIC     *cyclic;

528:   PetscNewLog(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->reset                = SVDReset_Cyclic;
535:   svd->ops->view                 = SVDView_Cyclic;
536:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
537:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
538:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
539:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
540:   return(0);
541: }