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-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:   PetscBool explicitmatrix;
 34:   EPS       eps;
 35:   PetscBool setfromoptionscalled;
 36:   Mat       mat;
 37:   Vec       x1,x2,y1,y2;
 38: } SVD_CYCLIC;

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

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

 74: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
 75: {

 79:   VecSet(diag,0.0);
 80:   return(0);
 81: }

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

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

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

189:   if (svd->ncv != svd->n) {
190:     VecDestroyVecs(svd->n,&svd->U);
191:     VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
192:     PetscLogObjectParents(svd,svd->ncv,svd->U);
193:   }
194:   return(0);
195: }

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

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

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

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

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

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

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

280:   cyclic->setfromoptionscalled = PETSC_TRUE;
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(svd,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(svd,cyclic->eps);
434:     if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
435:     EPSSetIP(cyclic->eps,svd->ip);
436:     EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
437:     EPSMonitorSet(cyclic->eps,SVDMonitor_Cyclic,svd,NULL);
438:   }
439:   *eps = cyclic->eps;
440:   return(0);
441: }

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

449:    Not Collective

451:    Input Parameter:
452: .  svd - singular value solver

454:    Output Parameter:
455: .  eps - the eigensolver object

457:    Level: advanced

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

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

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

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

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

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

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

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

524: PETSC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
525: {
527:   SVD_CYCLIC     *cyclic;

530:   PetscNewLog(svd,SVD_CYCLIC,&cyclic);
531:   svd->data                      = (void*)cyclic;
532:   svd->ops->solve                = SVDSolve_Cyclic;
533:   svd->ops->setup                = SVDSetUp_Cyclic;
534:   svd->ops->setfromoptions       = SVDSetFromOptions_Cyclic;
535:   svd->ops->destroy              = SVDDestroy_Cyclic;
536:   svd->ops->reset                = SVDReset_Cyclic;
537:   svd->ops->view                 = SVDView_Cyclic;
538:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
539:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
540:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
541:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
542:   return(0);
543: }