Actual source code: cyclic.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc singular value solver: "cyclic"

 13:    Method: Uses a Hermitian eigensolver for H(A) = [ 0  A ; A^T 0 ]
 14: */

 16: #include <slepc/private/svdimpl.h>
 17: #include <slepc/private/epsimpl.h>
 18: #include "cyclic.h"

 20: static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
 21: {
 22:   SVD_CYCLIC_SHELL  *ctx;
 23:   const PetscScalar *px;
 24:   PetscScalar       *py;
 25:   PetscInt          m;

 27:   MatShellGetContext(B,&ctx);
 28:   MatGetLocalSize(ctx->A,&m,NULL);
 29:   VecGetArrayRead(x,&px);
 30:   VecGetArrayWrite(y,&py);
 31:   VecPlaceArray(ctx->x1,px);
 32:   VecPlaceArray(ctx->x2,px+m);
 33:   VecPlaceArray(ctx->y1,py);
 34:   VecPlaceArray(ctx->y2,py+m);
 35:   MatMult(ctx->A,ctx->x2,ctx->y1);
 36:   MatMult(ctx->AT,ctx->x1,ctx->y2);
 37:   VecResetArray(ctx->x1);
 38:   VecResetArray(ctx->x2);
 39:   VecResetArray(ctx->y1);
 40:   VecResetArray(ctx->y2);
 41:   VecRestoreArrayRead(x,&px);
 42:   VecRestoreArrayWrite(y,&py);
 43:   PetscFunctionReturn(0);
 44: }

 46: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
 47: {
 48:   VecSet(diag,0.0);
 49:   PetscFunctionReturn(0);
 50: }

 52: static PetscErrorCode MatDestroy_Cyclic(Mat B)
 53: {
 54:   SVD_CYCLIC_SHELL *ctx;

 56:   MatShellGetContext(B,&ctx);
 57:   VecDestroy(&ctx->x1);
 58:   VecDestroy(&ctx->x2);
 59:   VecDestroy(&ctx->y1);
 60:   VecDestroy(&ctx->y2);
 61:   PetscFree(ctx);
 62:   PetscFunctionReturn(0);
 63: }

 65: /*
 66:    Builds cyclic matrix   C = | 0   A |
 67:                               | AT  0 |
 68: */
 69: static PetscErrorCode SVDCyclicGetCyclicMat(SVD svd,Mat A,Mat AT,Mat *C)
 70: {
 71:   SVD_CYCLIC       *cyclic = (SVD_CYCLIC*)svd->data;
 72:   SVD_CYCLIC_SHELL *ctx;
 73:   PetscInt         i,M,N,m,n,Istart,Iend;
 74:   VecType          vtype;
 75:   Mat              Zm,Zn;
 76: #if defined(PETSC_HAVE_CUDA)
 77:   PetscBool        cuda;
 78: #endif

 80:   MatGetSize(A,&M,&N);
 81:   MatGetLocalSize(A,&m,&n);

 83:   if (cyclic->explicitmatrix) {
 85:     MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
 86:     MatSetSizes(Zm,m,m,M,M);
 87:     MatSetFromOptions(Zm);
 88:     MatSetUp(Zm);
 89:     MatGetOwnershipRange(Zm,&Istart,&Iend);
 90:     for (i=Istart;i<Iend;i++) MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
 91:     MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
 92:     MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
 93:     MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
 94:     MatSetSizes(Zn,n,n,N,N);
 95:     MatSetFromOptions(Zn);
 96:     MatSetUp(Zn);
 97:     MatGetOwnershipRange(Zn,&Istart,&Iend);
 98:     for (i=Istart;i<Iend;i++) MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
 99:     MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
100:     MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
101:     MatCreateTile(1.0,Zm,1.0,A,1.0,AT,1.0,Zn,C);
102:     MatDestroy(&Zm);
103:     MatDestroy(&Zn);
104:   } else {
105:     PetscNew(&ctx);
106:     ctx->A       = A;
107:     ctx->AT      = AT;
108:     ctx->swapped = svd->swapped;
109:     MatCreateVecsEmpty(A,&ctx->x2,&ctx->x1);
110:     MatCreateVecsEmpty(A,&ctx->y2,&ctx->y1);
111:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x1);
112:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x2);
113:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y1);
114:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y2);
115:     MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C);
116:     MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic);
117:     MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cyclic);
118: #if defined(PETSC_HAVE_CUDA)
119:     PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&cuda,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
120:     if (cuda) MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic_CUDA);
121:     else
122: #endif
123:       MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic);
124:     MatGetVecType(A,&vtype);
125:     MatSetVecType(*C,vtype);
126:   }
127:   PetscLogObjectParent((PetscObject)svd,(PetscObject)*C);
128:   PetscFunctionReturn(0);
129: }

131: static PetscErrorCode MatMult_ECross(Mat B,Vec x,Vec y)
132: {
133:   SVD_CYCLIC_SHELL  *ctx;
134:   const PetscScalar *px;
135:   PetscScalar       *py;
136:   PetscInt          mn,m,n;

138:   MatShellGetContext(B,&ctx);
139:   MatGetLocalSize(ctx->A,NULL,&n);
140:   VecGetLocalSize(y,&mn);
141:   m = mn-n;
142:   VecGetArrayRead(x,&px);
143:   VecGetArrayWrite(y,&py);
144:   VecPlaceArray(ctx->x1,px);
145:   VecPlaceArray(ctx->x2,px+m);
146:   VecPlaceArray(ctx->y1,py);
147:   VecPlaceArray(ctx->y2,py+m);
148:   VecCopy(ctx->x1,ctx->y1);
149:   MatMult(ctx->A,ctx->x2,ctx->w);
150:   MatMult(ctx->AT,ctx->w,ctx->y2);
151:   VecResetArray(ctx->x1);
152:   VecResetArray(ctx->x2);
153:   VecResetArray(ctx->y1);
154:   VecResetArray(ctx->y2);
155:   VecRestoreArrayRead(x,&px);
156:   VecRestoreArrayWrite(y,&py);
157:   PetscFunctionReturn(0);
158: }

160: static PetscErrorCode MatGetDiagonal_ECross(Mat B,Vec d)
161: {
162:   SVD_CYCLIC_SHELL  *ctx;
163:   PetscScalar       *pd;
164:   PetscMPIInt       len;
165:   PetscInt          mn,m,n,N,i,j,start,end,ncols;
166:   PetscScalar       *work1,*work2,*diag;
167:   const PetscInt    *cols;
168:   const PetscScalar *vals;

170:   MatShellGetContext(B,&ctx);
171:   MatGetLocalSize(ctx->A,NULL,&n);
172:   VecGetLocalSize(d,&mn);
173:   m = mn-n;
174:   VecGetArrayWrite(d,&pd);
175:   VecPlaceArray(ctx->y1,pd);
176:   VecSet(ctx->y1,1.0);
177:   VecResetArray(ctx->y1);
178:   VecPlaceArray(ctx->y2,pd+m);
179:   if (!ctx->diag) {
180:     /* compute diagonal from rows and store in ctx->diag */
181:     VecDuplicate(ctx->y2,&ctx->diag);
182:     MatGetSize(ctx->A,NULL,&N);
183:     PetscCalloc2(N,&work1,N,&work2);
184:     if (ctx->swapped) {
185:       MatGetOwnershipRange(ctx->AT,&start,&end);
186:       for (i=start;i<end;i++) {
187:         MatGetRow(ctx->AT,i,&ncols,NULL,&vals);
188:         for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
189:         MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals);
190:       }
191:     } else {
192:       MatGetOwnershipRange(ctx->A,&start,&end);
193:       for (i=start;i<end;i++) {
194:         MatGetRow(ctx->A,i,&ncols,&cols,&vals);
195:         for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
196:         MatRestoreRow(ctx->A,i,&ncols,&cols,&vals);
197:       }
198:     }
199:     PetscMPIIntCast(N,&len);
200:     MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B));
201:     VecGetOwnershipRange(ctx->diag,&start,&end);
202:     VecGetArrayWrite(ctx->diag,&diag);
203:     for (i=start;i<end;i++) diag[i-start] = work2[i];
204:     VecRestoreArrayWrite(ctx->diag,&diag);
205:     PetscFree2(work1,work2);
206:   }
207:   VecCopy(ctx->diag,ctx->y2);
208:   VecResetArray(ctx->y2);
209:   VecRestoreArrayWrite(d,&pd);
210:   PetscFunctionReturn(0);
211: }

213: static PetscErrorCode MatDestroy_ECross(Mat B)
214: {
215:   SVD_CYCLIC_SHELL *ctx;

217:   MatShellGetContext(B,&ctx);
218:   VecDestroy(&ctx->x1);
219:   VecDestroy(&ctx->x2);
220:   VecDestroy(&ctx->y1);
221:   VecDestroy(&ctx->y2);
222:   VecDestroy(&ctx->diag);
223:   VecDestroy(&ctx->w);
224:   PetscFree(ctx);
225:   PetscFunctionReturn(0);
226: }

228: /*
229:    Builds extended cross product matrix   C = | I_m   0  |
230:                                               |  0  AT*A |
231:    t is an auxiliary Vec used to take the dimensions of the upper block
232: */
233: static PetscErrorCode SVDCyclicGetECrossMat(SVD svd,Mat A,Mat AT,Mat *C,Vec t)
234: {
235:   SVD_CYCLIC       *cyclic = (SVD_CYCLIC*)svd->data;
236:   SVD_CYCLIC_SHELL *ctx;
237:   PetscInt         i,M,N,m,n,Istart,Iend;
238:   VecType          vtype;
239:   Mat              Id,Zm,Zn,ATA;
240: #if defined(PETSC_HAVE_CUDA)
241:   PetscBool        cuda;
242: #endif

244:   MatGetSize(A,NULL,&N);
245:   MatGetLocalSize(A,NULL,&n);
246:   VecGetSize(t,&M);
247:   VecGetLocalSize(t,&m);

249:   if (cyclic->explicitmatrix) {
251:     MatCreateConstantDiagonal(PetscObjectComm((PetscObject)svd),m,m,M,M,1.0,&Id);
252:     MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
253:     MatSetSizes(Zm,m,n,M,N);
254:     MatSetFromOptions(Zm);
255:     MatSetUp(Zm);
256:     MatGetOwnershipRange(Zm,&Istart,&Iend);
257:     for (i=Istart;i<Iend;i++) {
258:       if (i<N) MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
259:     }
260:     MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
261:     MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
262:     MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
263:     MatSetSizes(Zn,n,m,N,M);
264:     MatSetFromOptions(Zn);
265:     MatSetUp(Zn);
266:     MatGetOwnershipRange(Zn,&Istart,&Iend);
267:     for (i=Istart;i<Iend;i++) {
268:       if (i<m) MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
269:     }
270:     MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
271:     MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
272:     MatProductCreate(AT,A,NULL,&ATA);
273:     MatProductSetType(ATA,MATPRODUCT_AB);
274:     MatProductSetFromOptions(ATA);
275:     MatProductSymbolic(ATA);
276:     MatProductNumeric(ATA);
277:     MatCreateTile(1.0,Id,1.0,Zm,1.0,Zn,1.0,ATA,C);
278:     MatDestroy(&Id);
279:     MatDestroy(&Zm);
280:     MatDestroy(&Zn);
281:     MatDestroy(&ATA);
282:   } else {
283:     PetscNew(&ctx);
284:     ctx->A       = A;
285:     ctx->AT      = AT;
286:     ctx->swapped = svd->swapped;
287:     VecDuplicateEmpty(t,&ctx->x1);
288:     VecDuplicateEmpty(t,&ctx->y1);
289:     MatCreateVecsEmpty(A,&ctx->x2,NULL);
290:     MatCreateVecsEmpty(A,&ctx->y2,NULL);
291:     MatCreateVecs(A,NULL,&ctx->w);
292:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x1);
293:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x2);
294:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y1);
295:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y2);
296:     MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C);
297:     MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ECross);
298:     MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_ECross);
299: #if defined(PETSC_HAVE_CUDA)
300:     PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&cuda,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
301:     if (cuda) MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross_CUDA);
302:     else
303: #endif
304:       MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross);
305:     MatGetVecType(A,&vtype);
306:     MatSetVecType(*C,vtype);
307:   }
308:   PetscLogObjectParent((PetscObject)svd,(PetscObject)*C);
309:   PetscFunctionReturn(0);
310: }

312: /* Convergence test relative to the norm of R (used in GSVD only) */
313: static PetscErrorCode EPSConv_Cyclic(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
314: {
315:   SVD svd = (SVD)ctx;

317:   *errest = res/PetscMax(svd->nrma,svd->nrmb);
318:   PetscFunctionReturn(0);
319: }

321: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
322: {
323:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC*)svd->data;
324:   PetscInt          M,N,m,n,p,k,i,isl,offset;
325:   const PetscScalar *isa;
326:   PetscScalar       *va;
327:   PetscBool         trackall,issinv;
328:   Vec               v,t;
329:   ST                st;

331:   MatGetSize(svd->A,&M,&N);
332:   MatGetLocalSize(svd->A,&m,&n);
333:   if (!cyclic->eps) SVDCyclicGetEPS(svd,&cyclic->eps);
334:   MatDestroy(&cyclic->C);
335:   MatDestroy(&cyclic->D);
336:   if (svd->isgeneralized) {
337:     if (svd->which==SVD_SMALLEST) {  /* alternative pencil */
338:       MatCreateVecs(svd->B,NULL,&t);
339:       SVDCyclicGetCyclicMat(svd,svd->B,svd->BT,&cyclic->C);
340:       SVDCyclicGetECrossMat(svd,svd->A,svd->AT,&cyclic->D,t);
341:     } else {
342:       MatCreateVecs(svd->A,NULL,&t);
343:       SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C);
344:       SVDCyclicGetECrossMat(svd,svd->B,svd->BT,&cyclic->D,t);
345:     }
346:     VecDestroy(&t);
347:     EPSSetOperators(cyclic->eps,cyclic->C,cyclic->D);
348:     EPSSetProblemType(cyclic->eps,EPS_GHEP);
349:   } else {
350:     SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C);
351:     EPSSetOperators(cyclic->eps,cyclic->C,NULL);
352:     EPSSetProblemType(cyclic->eps,EPS_HEP);
353:   }
354:   if (!cyclic->usereps) {
355:     if (svd->which == SVD_LARGEST) {
356:       EPSGetST(cyclic->eps,&st);
357:       PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv);
358:       if (issinv) EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE);
359:       else EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
360:     } else {
361:       if (svd->isgeneralized) {  /* computes sigma^{-1} via alternative pencil */
362:         EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
363:       } else {
364:         EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
365:         EPSSetTarget(cyclic->eps,0.0);
366:       }
367:     }
368:     EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
369:     EPSSetTolerances(cyclic->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it);
370:     switch (svd->conv) {
371:     case SVD_CONV_ABS:
372:       EPSSetConvergenceTest(cyclic->eps,EPS_CONV_ABS);break;
373:     case SVD_CONV_REL:
374:       EPSSetConvergenceTest(cyclic->eps,EPS_CONV_REL);break;
375:     case SVD_CONV_NORM:
376:       if (svd->isgeneralized) {
377:         if (!svd->nrma) MatNorm(svd->OP,NORM_INFINITY,&svd->nrma);
378:         if (!svd->nrmb) MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb);
379:         EPSSetConvergenceTestFunction(cyclic->eps,EPSConv_Cyclic,svd,NULL);
380:       } else {
381:         EPSSetConvergenceTest(cyclic->eps,EPS_CONV_NORM);break;
382:       }
383:       break;
384:     case SVD_CONV_MAXIT:
385:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
386:     case SVD_CONV_USER:
387:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
388:     }
389:   }
390:   SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
391:   /* Transfer the trackall option from svd to eps */
392:   SVDGetTrackAll(svd,&trackall);
393:   EPSSetTrackAll(cyclic->eps,trackall);
394:   /* Transfer the initial subspace from svd to eps */
395:   if (svd->nini<0 || svd->ninil<0) {
396:     for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
397:       MatCreateVecs(cyclic->C,&v,NULL);
398:       VecGetArrayWrite(v,&va);
399:       if (svd->isgeneralized) MatGetLocalSize(svd->B,&p,NULL);
400:       k = (svd->isgeneralized && svd->which==SVD_SMALLEST)? p: m;  /* size of upper block row */
401:       if (i<-svd->ninil) {
402:         VecGetArrayRead(svd->ISL[i],&isa);
403:         if (svd->isgeneralized) {
404:           VecGetLocalSize(svd->ISL[i],&isl);
406:           offset = (svd->which==SVD_SMALLEST)? m: 0;
407:           PetscArraycpy(va,isa+offset,k);
408:         } else {
409:           VecGetLocalSize(svd->ISL[i],&isl);
411:           PetscArraycpy(va,isa,k);
412:         }
413:         VecRestoreArrayRead(svd->IS[i],&isa);
414:       } else PetscArrayzero(&va,k);
415:       if (i<-svd->nini) {
416:         VecGetLocalSize(svd->IS[i],&isl);
418:         VecGetArrayRead(svd->IS[i],&isa);
419:         PetscArraycpy(va+k,isa,n);
420:         VecRestoreArrayRead(svd->IS[i],&isa);
421:       } else PetscArrayzero(va+k,n);
422:       VecRestoreArrayWrite(v,&va);
423:       VecDestroy(&svd->IS[i]);
424:       svd->IS[i] = v;
425:     }
426:     svd->nini = PetscMin(svd->nini,svd->ninil);
427:     EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
428:     SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
429:     SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
430:   }
431:   EPSSetUp(cyclic->eps);
432:   EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
433:   svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
434:   EPSGetTolerances(cyclic->eps,NULL,&svd->max_it);
435:   if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;

437:   svd->leftbasis = PETSC_TRUE;
438:   SVDAllocateSolution(svd,0);
439:   PetscFunctionReturn(0);
440: }

442: PetscErrorCode SVDSolve_Cyclic(SVD svd)
443: {
444:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;
445:   PetscInt       i,j,nconv;
446:   PetscScalar    sigma;

448:   EPSSolve(cyclic->eps);
449:   EPSGetConverged(cyclic->eps,&nconv);
450:   EPSGetIterationNumber(cyclic->eps,&svd->its);
451:   EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
452:   for (i=0,j=0;i<nconv;i++) {
453:     EPSGetEigenvalue(cyclic->eps,i,&sigma,NULL);
454:     if (PetscRealPart(sigma) > 0.0) {
455:       if (svd->isgeneralized && svd->which==SVD_SMALLEST) svd->sigma[j] = 1.0/PetscRealPart(sigma);
456:       else svd->sigma[j] = PetscRealPart(sigma);
457:       j++;
458:     }
459:   }
460:   svd->nconv = j;
461:   PetscFunctionReturn(0);
462: }

464: PetscErrorCode SVDComputeVectors_Cyclic(SVD svd)
465: {
466:   SVD_CYCLIC        *cyclic = (SVD_CYCLIC*)svd->data;
467:   PetscInt          i,j,m,p,nconv;
468:   PetscScalar       *dst,sigma;
469:   const PetscScalar *src,*px;
470:   Vec               u,v,x,x1,x2,uv;

472:   EPSGetConverged(cyclic->eps,&nconv);
473:   MatCreateVecs(cyclic->C,&x,NULL);
474:   MatGetLocalSize(svd->A,&m,NULL);
475:   if (svd->isgeneralized && svd->which==SVD_SMALLEST) MatCreateVecsEmpty(svd->B,&x1,&x2);
476:   else MatCreateVecsEmpty(svd->A,&x2,&x1);
477:   if (svd->isgeneralized) {
478:     MatCreateVecs(svd->A,NULL,&u);
479:     MatCreateVecs(svd->B,NULL,&v);
480:     MatGetLocalSize(svd->B,&p,NULL);
481:   }
482:   for (i=0,j=0;i<nconv;i++) {
483:     EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
484:     if (PetscRealPart(sigma) > 0.0) {
485:       if (svd->isgeneralized) {
486:         if (svd->which==SVD_SMALLEST) {
487:           /* evec_i = 1/sqrt(2)*[ v_i; w_i ],  w_i = x_i/c_i */
488:           VecGetArrayRead(x,&px);
489:           VecPlaceArray(x2,px);
490:           VecPlaceArray(x1,px+p);
491:           VecCopy(x2,v);
492:           VecScale(v,PETSC_SQRT2);  /* v_i = sqrt(2)*evec_i_1 */
493:           VecScale(x1,PETSC_SQRT2); /* w_i = sqrt(2)*evec_i_2 */
494:           MatMult(svd->A,x1,u);     /* A*w_i = u_i */
495:           VecScale(x1,1.0/PetscSqrtScalar(1.0+sigma*sigma));  /* x_i = w_i*c_i */
496:           BVInsertVec(svd->V,j,x1);
497:           VecResetArray(x2);
498:           VecResetArray(x1);
499:           VecRestoreArrayRead(x,&px);
500:         } else {
501:           /* evec_i = 1/sqrt(2)*[ u_i; w_i ],  w_i = x_i/s_i */
502:           VecGetArrayRead(x,&px);
503:           VecPlaceArray(x1,px);
504:           VecPlaceArray(x2,px+m);
505:           VecCopy(x1,u);
506:           VecScale(u,PETSC_SQRT2);  /* u_i = sqrt(2)*evec_i_1 */
507:           VecScale(x2,PETSC_SQRT2); /* w_i = sqrt(2)*evec_i_2 */
508:           MatMult(svd->B,x2,v);     /* B*w_i = v_i */
509:           VecScale(x2,1.0/PetscSqrtScalar(1.0+sigma*sigma));  /* x_i = w_i*s_i */
510:           BVInsertVec(svd->V,j,x2);
511:           VecResetArray(x1);
512:           VecResetArray(x2);
513:           VecRestoreArrayRead(x,&px);
514:         }
515:         /* copy [u;v] to U[j] */
516:         BVGetColumn(svd->U,j,&uv);
517:         VecGetArrayWrite(uv,&dst);
518:         VecGetArrayRead(u,&src);
519:         PetscArraycpy(dst,src,m);
520:         VecRestoreArrayRead(u,&src);
521:         VecGetArrayRead(v,&src);
522:         PetscArraycpy(dst+m,src,p);
523:         VecRestoreArrayRead(v,&src);
524:         VecRestoreArrayWrite(uv,&dst);
525:         BVRestoreColumn(svd->U,j,&uv);
526:       } else {
527:         VecGetArrayRead(x,&px);
528:         VecPlaceArray(x1,px);
529:         VecPlaceArray(x2,px+m);
530:         BVInsertVec(svd->U,j,x1);
531:         BVScaleColumn(svd->U,j,PETSC_SQRT2);
532:         BVInsertVec(svd->V,j,x2);
533:         BVScaleColumn(svd->V,j,PETSC_SQRT2);
534:         VecResetArray(x1);
535:         VecResetArray(x2);
536:         VecRestoreArrayRead(x,&px);
537:       }
538:       j++;
539:     }
540:   }
541:   VecDestroy(&x);
542:   VecDestroy(&x1);
543:   VecDestroy(&x2);
544:   if (svd->isgeneralized) {
545:     VecDestroy(&u);
546:     VecDestroy(&v);
547:   }
548:   PetscFunctionReturn(0);
549: }

551: static PetscErrorCode EPSMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
552: {
553:   PetscInt       i,j;
554:   SVD            svd = (SVD)ctx;
555:   PetscScalar    er,ei;

557:   nconv = 0;
558:   for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
559:     er = eigr[i]; ei = eigi[i];
560:     STBackTransform(eps->st,1,&er,&ei);
561:     if (PetscRealPart(er) > 0.0) {
562:       svd->sigma[j] = PetscRealPart(er);
563:       svd->errest[j] = errest[i];
564:       if (errest[i] && errest[i] < svd->tol) nconv++;
565:       j++;
566:     }
567:   }
568:   nest = j;
569:   SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
570:   PetscFunctionReturn(0);
571: }

573: PetscErrorCode SVDSetFromOptions_Cyclic(PetscOptionItems *PetscOptionsObject,SVD svd)
574: {
575:   PetscBool      set,val;
576:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;
577:   ST             st;

579:   PetscOptionsHead(PetscOptionsObject,"SVD Cyclic Options");

581:     PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
582:     if (set) SVDCyclicSetExplicitMatrix(svd,val);

584:   PetscOptionsTail();

586:   if (!cyclic->eps) SVDCyclicGetEPS(svd,&cyclic->eps);
587:   if (!cyclic->explicitmatrix && !cyclic->usereps) {
588:     /* use as default an ST with shell matrix and Jacobi */
589:     EPSGetST(cyclic->eps,&st);
590:     STSetMatMode(st,ST_MATMODE_SHELL);
591:   }
592:   EPSSetFromOptions(cyclic->eps);
593:   PetscFunctionReturn(0);
594: }

596: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmat)
597: {
598:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

600:   if (cyclic->explicitmatrix != explicitmat) {
601:     cyclic->explicitmatrix = explicitmat;
602:     svd->state = SVD_STATE_INITIAL;
603:   }
604:   PetscFunctionReturn(0);
605: }

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

611:    Logically Collective on svd

613:    Input Parameters:
614: +  svd         - singular value solver
615: -  explicitmat - boolean flag indicating if H(A) is built explicitly

617:    Options Database Key:
618: .  -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag

620:    Level: advanced

622: .seealso: SVDCyclicGetExplicitMatrix()
623: @*/
624: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmat)
625: {
628:   PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
629:   PetscFunctionReturn(0);
630: }

632: static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmat)
633: {
634:   SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;

636:   *explicitmat = cyclic->explicitmatrix;
637:   PetscFunctionReturn(0);
638: }

640: /*@
641:    SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly.

643:    Not Collective

645:    Input Parameter:
646: .  svd  - singular value solver

648:    Output Parameter:
649: .  explicitmat - the mode flag

651:    Level: advanced

653: .seealso: SVDCyclicSetExplicitMatrix()
654: @*/
655: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
656: {
659:   PetscUseMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
660:   PetscFunctionReturn(0);
661: }

663: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
664: {
665:   SVD_CYCLIC      *cyclic = (SVD_CYCLIC*)svd->data;

667:   PetscObjectReference((PetscObject)eps);
668:   EPSDestroy(&cyclic->eps);
669:   cyclic->eps = eps;
670:   cyclic->usereps = PETSC_TRUE;
671:   PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
672:   svd->state = SVD_STATE_INITIAL;
673:   PetscFunctionReturn(0);
674: }

676: /*@
677:    SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
678:    singular value solver.

680:    Collective on svd

682:    Input Parameters:
683: +  svd - singular value solver
684: -  eps - the eigensolver object

686:    Level: advanced

688: .seealso: SVDCyclicGetEPS()
689: @*/
690: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
691: {
695:   PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
696:   PetscFunctionReturn(0);
697: }

699: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
700: {
701:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

703:   if (!cyclic->eps) {
704:     EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
705:     PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
706:     EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
707:     EPSAppendOptionsPrefix(cyclic->eps,"svd_cyclic_");
708:     PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
709:     PetscObjectSetOptions((PetscObject)cyclic->eps,((PetscObject)svd)->options);
710:     EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
711:     EPSMonitorSet(cyclic->eps,EPSMonitor_Cyclic,svd,NULL);
712:   }
713:   *eps = cyclic->eps;
714:   PetscFunctionReturn(0);
715: }

717: /*@
718:    SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
719:    to the singular value solver.

721:    Not Collective

723:    Input Parameter:
724: .  svd - singular value solver

726:    Output Parameter:
727: .  eps - the eigensolver object

729:    Level: advanced

731: .seealso: SVDCyclicSetEPS()
732: @*/
733: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
734: {
737:   PetscUseMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
738:   PetscFunctionReturn(0);
739: }

741: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
742: {
743:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;
744:   PetscBool      isascii;

746:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
747:   if (isascii) {
748:     if (!cyclic->eps) SVDCyclicGetEPS(svd,&cyclic->eps);
749:     PetscViewerASCIIPrintf(viewer,"  %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
750:     PetscViewerASCIIPushTab(viewer);
751:     EPSView(cyclic->eps,viewer);
752:     PetscViewerASCIIPopTab(viewer);
753:   }
754:   PetscFunctionReturn(0);
755: }

757: PetscErrorCode SVDReset_Cyclic(SVD svd)
758: {
759:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

761:   EPSReset(cyclic->eps);
762:   MatDestroy(&cyclic->C);
763:   MatDestroy(&cyclic->D);
764:   PetscFunctionReturn(0);
765: }

767: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
768: {
769:   SVD_CYCLIC     *cyclic = (SVD_CYCLIC*)svd->data;

771:   EPSDestroy(&cyclic->eps);
772:   PetscFree(svd->data);
773:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
774:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
775:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
776:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
777:   PetscFunctionReturn(0);
778: }

780: SLEPC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
781: {
782:   SVD_CYCLIC     *cyclic;

784:   PetscNewLog(svd,&cyclic);
785:   svd->data                = (void*)cyclic;
786:   svd->ops->solve          = SVDSolve_Cyclic;
787:   svd->ops->solveg         = SVDSolve_Cyclic;
788:   svd->ops->setup          = SVDSetUp_Cyclic;
789:   svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
790:   svd->ops->destroy        = SVDDestroy_Cyclic;
791:   svd->ops->reset          = SVDReset_Cyclic;
792:   svd->ops->view           = SVDView_Cyclic;
793:   svd->ops->computevectors = SVDComputeVectors_Cyclic;
794:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
795:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
796:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
797:   PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
798:   PetscFunctionReturn(0);
799: }