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-2011, Universitat 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>                /*I "slepcsvd.h" I*/
 30: #include <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 ShellMatMult_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;
 50: 
 52:   MatShellGetContext(B,(void**)&svd);
 53:   cyclic = (SVD_CYCLIC*)svd->data;
 54:   SVDMatGetLocalSize(svd,&m,PETSC_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 ShellMatGetDiagonal_Cyclic(Mat B,Vec diag)
 75: {
 77: 
 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(((PetscObject)svd)->comm,PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
102:       MatCreate(((PetscObject)svd)->comm,&Zm);
103:       MatSetSizes(Zm,m,m,M,M);
104:       MatSetFromOptions(Zm);
105:       MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
106:       MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
107:       MatCreate(((PetscObject)svd)->comm,&Zn);
108:       MatSetSizes(Zn,n,n,N,N);
109:       MatSetFromOptions(Zn);
110:       MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
111:       MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
112:       SlepcMatTile(0.0,Zm,1.0,svd->A,1.0,svd->AT,0.0,Zn,&cyclic->mat);
113:       MatDestroy(&Zm);
114:       MatDestroy(&Zn);
115:     } else {
116:       VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->x1);
117:       VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->x2);
118:       VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&cyclic->y1);
119:       VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&cyclic->y2);
120:       MatCreateShell(((PetscObject)svd)->comm,m+n,m+n,M+N,M+N,svd,&cyclic->mat);
121:       MatShellSetOperation(cyclic->mat,MATOP_MULT,(void(*)(void))ShellMatMult_Cyclic);
122:       MatShellSetOperation(cyclic->mat,MATOP_GET_DIAGONAL,(void(*)(void))ShellMatGetDiagonal_Cyclic);
123:     }
124:   }

126:   EPSSetOperators(cyclic->eps,cyclic->mat,PETSC_NULL);
127:   EPSSetProblemType(cyclic->eps,EPS_HEP);
128:   EPSSetWhichEigenpairs(cyclic->eps,svd->which == SVD_LARGEST ? EPS_LARGEST_REAL : EPS_SMALLEST_MAGNITUDE);
129:   EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
130:   EPSSetTolerances(cyclic->eps,svd->tol,svd->max_it);
131:   /* Transfer the trackall option from svd to eps */
132:   SVDGetTrackAll(svd,&trackall);
133:   EPSSetTrackAll(cyclic->eps,trackall);
134:   /* Transfer the initial subspace from svd to eps */
135:   if (svd->nini < 0) {
136:     for (i=0; i<-svd->nini; i++) {
137:       MatGetVecs(cyclic->mat,&v,PETSC_NULL);
138:       VecGetArray(v,&va);
139:       VecGetArrayRead(svd->IS[i],&isa);
140:       VecGetSize(svd->IS[i],&isl);
141:       if (isl == m) {
142:         PetscMemcpy(va,isa,sizeof(PetscScalar)*m);
143:         PetscMemzero(&va[m],sizeof(PetscScalar)*n);
144:       } else if (isl == n) {
145:         PetscMemzero(va,sizeof(PetscScalar)*m);
146:         PetscMemcpy(&va[m],isa,sizeof(PetscScalar)*n);
147:       } else {
148:         SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_SUP,"Size of the initial subspace vectors should match to some dimension of A");
149:       }
150:       VecRestoreArray(v,&va);
151:       VecRestoreArrayRead(svd->IS[i],&isa);
152:       VecDestroy(&svd->IS[i]);
153:       svd->IS[i] = v;
154:     }
155:     EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
156:     for (i=0; i<-svd->nini; i++) {
157:       VecDestroy(&svd->IS[i]);
158:     }
159:     PetscFree(svd->IS);
160:     svd->nini = 0;
161:   }
162:   if (cyclic->setfromoptionscalled) {
163:     EPSSetFromOptions(cyclic->eps);
164:     cyclic->setfromoptionscalled = PETSC_FALSE;
165:   }
166:   EPSSetUp(cyclic->eps);
167:   EPSGetDimensions(cyclic->eps,PETSC_NULL,&svd->ncv,&svd->mpd);
168:   svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
169:   EPSGetTolerances(cyclic->eps,&svd->tol,&svd->max_it);

171:   if (svd->ncv != svd->n) {
172:     VecDestroyVecs(svd->n,&svd->U);
173:     VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
174:   }
175:   return(0);
176: }

180: PetscErrorCode SVDSolve_Cyclic(SVD svd)
181: {
182:   PetscErrorCode    ierr;
183:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC*)svd->data;
184:   PetscInt          i,j,M,N,m,n;
185:   PetscScalar       sigma;
186:   const PetscScalar *px;
187:   Vec               x,x1,x2;
188: 
190:   EPSSolve(cyclic->eps);
191:   EPSGetConverged(cyclic->eps,&svd->nconv);
192:   EPSGetIterationNumber(cyclic->eps,&svd->its);
193:   EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);

195:   MatGetVecs(cyclic->mat,&x,PETSC_NULL);
196:   SVDMatGetSize(svd,&M,&N);
197:   SVDMatGetLocalSize(svd,&m,&n);
198:   VecCreateMPIWithArray(((PetscObject)svd)->comm,m,M,PETSC_NULL,&x1);
199:   VecCreateMPIWithArray(((PetscObject)svd)->comm,n,N,PETSC_NULL,&x2);
200:   for (i=0,j=0;i<svd->nconv;i++) {
201:     EPSGetEigenpair(cyclic->eps,i,&sigma,PETSC_NULL,x,PETSC_NULL);
202:     if (PetscRealPart(sigma) > 0.0) {
203:       svd->sigma[j] = PetscRealPart(sigma);
204:       VecGetArrayRead(x,&px);
205:       VecPlaceArray(x1,px);
206:       VecPlaceArray(x2,px+m);
207:       VecCopy(x1,svd->U[j]);
208:       VecScale(svd->U[j],1.0/PetscSqrtReal(2.0));
209:       VecCopy(x2,svd->V[j]);
210:       VecScale(svd->V[j],1.0/PetscSqrtReal(2.0));
211:       VecResetArray(x1);
212:       VecResetArray(x2);
213:       VecRestoreArrayRead(x,&px);
214:       j++;
215:     }
216:   }
217:   svd->nconv = j;

219:   VecDestroy(&x);
220:   VecDestroy(&x1);
221:   VecDestroy(&x2);
222:   return(0);
223: }

227: static PetscErrorCode SVDMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
228: {
229:   PetscInt       i,j;
230:   SVD            svd = (SVD)ctx;
231:   PetscScalar    er,ei;

235:   nconv = 0;
236:   for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
237:     er = eigr[i]; ei = eigi[i];
238:     STBackTransform(eps->OP,1,&er,&ei);
239:     if (PetscRealPart(er) > 0.0) {
240:       svd->sigma[j] = PetscRealPart(er);
241:       svd->errest[j] = errest[i];
242:       if (errest[i] < svd->tol) nconv++;
243:       j++;
244:     }
245:   }
246:   nest = j;
247:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
248:   return(0);
249: }

253: PetscErrorCode SVDSetFromOptions_Cyclic(SVD svd)
254: {
256:   PetscBool      set,val;
257:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;
258:   ST             st;

261:   cyclic->setfromoptionscalled = PETSC_TRUE;
262:   PetscOptionsHead("SVD Cyclic Options");
263:   PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
264:   if (set) {
265:     SVDCyclicSetExplicitMatrix(svd,val);
266:   }
267:   if (!cyclic->explicitmatrix) {
268:     /* use as default an ST with shell matrix and Jacobi */
269:     EPSGetST(cyclic->eps,&st);
270:     STSetMatMode(st,ST_MATMODE_SHELL);
271:   }
272:   PetscOptionsTail();
273:   return(0);
274: }

279: PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
280: {
281:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

284:   cyclic->explicitmatrix = explicitmatrix;
285:   return(0);
286: }

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

295:    Logically Collective on SVD

297:    Input Parameters:
298: +  svd      - singular value solver
299: -  explicit - boolean flag indicating if H(A) is built explicitly

301:    Options Database Key:
302: .  -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag

304:    Level: advanced

306: .seealso: SVDCyclicGetExplicitMatrix()
307: @*/
308: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
309: {

315:   PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
316:   return(0);
317: }

322: PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
323: {
324:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

327:   *explicitmatrix = cyclic->explicitmatrix;
328:   return(0);
329: }

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

337:    Not Collective

339:    Input Parameter:
340: .  svd  - singular value solver

342:    Output Parameter:
343: .  explicit - the mode flag

345:    Level: advanced

347: .seealso: SVDCyclicSetExplicitMatrix()
348: @*/
349: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
350: {

356:   PetscTryMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
357:   return(0);
358: }

363: PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
364: {
365:   PetscErrorCode  ierr;
366:   SVD_CYCLIC      *cyclic = (SVD_CYCLIC*)svd->data;

369:   PetscObjectReference((PetscObject)eps);
370:   EPSDestroy(&cyclic->eps);
371:   cyclic->eps = eps;
372:   PetscLogObjectParent(svd,cyclic->eps);
373:   svd->setupcalled = 0;
374:   return(0);
375: }

380: /*@
381:    SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
382:    singular value solver. 

384:    Collective on SVD

386:    Input Parameters:
387: +  svd - singular value solver
388: -  eps - the eigensolver object

390:    Level: advanced

392: .seealso: SVDCyclicGetEPS()
393: @*/
394: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
395: {

402:   PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
403:   return(0);
404: }

409: PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
410: {
411:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

414:   *eps = cyclic->eps;
415:   return(0);
416: }

421: /*@
422:    SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
423:    to the singular value solver.

425:    Not Collective

427:    Input Parameter:
428: .  svd - singular value solver

430:    Output Parameter:
431: .  eps - the eigensolver object

433:    Level: advanced

435: .seealso: SVDCyclicSetEPS()
436: @*/
437: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
438: {

444:   PetscTryMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
445:   return(0);
446: }

450: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
451: {
453:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

456:   PetscViewerASCIIPrintf(viewer,"  Cyclic: %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
457:   PetscViewerASCIIPushTab(viewer);
458:   EPSView(cyclic->eps,viewer);
459:   PetscViewerASCIIPopTab(viewer);
460:   return(0);
461: }

465: PetscErrorCode SVDReset_Cyclic(SVD svd)
466: {
468:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

471:   EPSReset(cyclic->eps);
472:   MatDestroy(&cyclic->mat);
473:   VecDestroy(&cyclic->x1);
474:   VecDestroy(&cyclic->x2);
475:   VecDestroy(&cyclic->y1);
476:   VecDestroy(&cyclic->y2);
477:   return(0);
478: }

482: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
483: {
485:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

488:   EPSDestroy(&cyclic->eps);
489:   PetscFree(svd->data);
490:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","",PETSC_NULL);
491:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","",PETSC_NULL);
492:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","",PETSC_NULL);
493:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","",PETSC_NULL);
494:   return(0);
495: }

500: PetscErrorCode SVDCreate_Cyclic(SVD svd)
501: {
503:   SVD_CYCLIC     *cyclic;
504: 
506:   PetscNewLog(svd,SVD_CYCLIC,&cyclic);
507:   svd->data                      = (void *)cyclic;
508:   svd->ops->solve                = SVDSolve_Cyclic;
509:   svd->ops->setup                = SVDSetUp_Cyclic;
510:   svd->ops->setfromoptions       = SVDSetFromOptions_Cyclic;
511:   svd->ops->destroy              = SVDDestroy_Cyclic;
512:   svd->ops->reset                = SVDReset_Cyclic;
513:   svd->ops->view                 = SVDView_Cyclic;
514:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetEPS_C","SVDCyclicSetEPS_Cyclic",SVDCyclicSetEPS_Cyclic);
515:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetEPS_C","SVDCyclicGetEPS_Cyclic",SVDCyclicGetEPS_Cyclic);
516:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C","SVDCyclicSetExplicitMatrix_Cyclic",SVDCyclicSetExplicitMatrix_Cyclic);
517:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C","SVDCyclicGetExplicitMatrix_Cyclic",SVDCyclicGetExplicitMatrix_Cyclic);

519:   EPSCreate(((PetscObject)svd)->comm,&cyclic->eps);
520:   EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
521:   EPSAppendOptionsPrefix(cyclic->eps,"svd_");
522:   PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
523:   PetscLogObjectParent(svd,cyclic->eps);
524:   if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
525:   EPSSetIP(cyclic->eps,svd->ip);
526:   EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
527:   EPSMonitorSet(cyclic->eps,SVDMonitor_Cyclic,svd,PETSC_NULL);
528:   cyclic->explicitmatrix = PETSC_FALSE;
529:   cyclic->mat = PETSC_NULL;
530:   cyclic->x1 = PETSC_NULL;
531:   cyclic->x2 = PETSC_NULL;
532:   cyclic->y1 = PETSC_NULL;
533:   cyclic->y2 = PETSC_NULL;
534:   cyclic->setfromoptionscalled = PETSC_FALSE;
535:   return(0);
536: }