Actual source code: cyclic.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: PetscErrorCode ierr;
23: SVD_CYCLIC_SHELL *ctx;
24: const PetscScalar *px;
25: PetscScalar *py;
26: PetscInt m;
29: MatShellGetContext(B,&ctx);
30: MatGetLocalSize(ctx->A,&m,NULL);
31: VecGetArrayRead(x,&px);
32: VecGetArrayWrite(y,&py);
33: VecPlaceArray(ctx->x1,px);
34: VecPlaceArray(ctx->x2,px+m);
35: VecPlaceArray(ctx->y1,py);
36: VecPlaceArray(ctx->y2,py+m);
37: MatMult(ctx->A,ctx->x2,ctx->y1);
38: MatMult(ctx->AT,ctx->x1,ctx->y2);
39: VecResetArray(ctx->x1);
40: VecResetArray(ctx->x2);
41: VecResetArray(ctx->y1);
42: VecResetArray(ctx->y2);
43: VecRestoreArrayRead(x,&px);
44: VecRestoreArrayWrite(y,&py);
45: return(0);
46: }
48: static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
49: {
53: VecSet(diag,0.0);
54: return(0);
55: }
57: static PetscErrorCode MatDestroy_Cyclic(Mat B)
58: {
59: PetscErrorCode ierr;
60: SVD_CYCLIC_SHELL *ctx;
63: MatShellGetContext(B,&ctx);
64: VecDestroy(&ctx->x1);
65: VecDestroy(&ctx->x2);
66: VecDestroy(&ctx->y1);
67: VecDestroy(&ctx->y2);
68: PetscFree(ctx);
69: return(0);
70: }
72: /*
73: Builds cyclic matrix C = | 0 A |
74: | AT 0 |
75: */
76: static PetscErrorCode SVDCyclicGetCyclicMat(SVD svd,Mat A,Mat AT,Mat *C)
77: {
78: PetscErrorCode ierr;
79: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
80: SVD_CYCLIC_SHELL *ctx;
81: PetscInt i,M,N,m,n,Istart,Iend;
82: VecType vtype;
83: Mat Zm,Zn;
84: #if defined(PETSC_HAVE_CUDA)
85: PetscBool cuda;
86: #endif
89: MatGetSize(A,&M,&N);
90: MatGetLocalSize(A,&m,&n);
92: if (cyclic->explicitmatrix) {
93: if (!svd->expltrans) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
94: MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
95: MatSetSizes(Zm,m,m,M,M);
96: MatSetFromOptions(Zm);
97: MatSetUp(Zm);
98: MatGetOwnershipRange(Zm,&Istart,&Iend);
99: for (i=Istart;i<Iend;i++) {
100: MatSetValue(Zm,i,i,0.0,INSERT_VALUES);
101: }
102: MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
104: MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
105: MatSetSizes(Zn,n,n,N,N);
106: MatSetFromOptions(Zn);
107: MatSetUp(Zn);
108: MatGetOwnershipRange(Zn,&Istart,&Iend);
109: for (i=Istart;i<Iend;i++) {
110: MatSetValue(Zn,i,i,0.0,INSERT_VALUES);
111: }
112: MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
114: MatCreateTile(1.0,Zm,1.0,A,1.0,AT,1.0,Zn,C);
115: MatDestroy(&Zm);
116: MatDestroy(&Zn);
117: } else {
118: PetscNew(&ctx);
119: ctx->A = A;
120: ctx->AT = AT;
121: ctx->swapped = svd->swapped;
122: MatCreateVecsEmpty(A,&ctx->x2,&ctx->x1);
123: MatCreateVecsEmpty(A,&ctx->y2,&ctx->y1);
124: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x1);
125: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x2);
126: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y1);
127: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y2);
128: MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C);
129: MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic);
130: MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cyclic);
131: #if defined(PETSC_HAVE_CUDA)
132: PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&cuda,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
133: if (cuda) {
134: MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic_CUDA);
135: } else
136: #endif
137: {
138: MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic);
139: }
140: MatGetVecType(A,&vtype);
141: MatSetVecType(*C,vtype);
142: }
143: PetscLogObjectParent((PetscObject)svd,(PetscObject)*C);
144: return(0);
145: }
147: static PetscErrorCode MatMult_ECross(Mat B,Vec x,Vec y)
148: {
149: PetscErrorCode ierr;
150: SVD_CYCLIC_SHELL *ctx;
151: const PetscScalar *px;
152: PetscScalar *py;
153: PetscInt mn,m,n;
156: MatShellGetContext(B,&ctx);
157: MatGetLocalSize(ctx->A,NULL,&n);
158: VecGetLocalSize(y,&mn);
159: m = mn-n;
160: VecGetArrayRead(x,&px);
161: VecGetArrayWrite(y,&py);
162: VecPlaceArray(ctx->x1,px);
163: VecPlaceArray(ctx->x2,px+m);
164: VecPlaceArray(ctx->y1,py);
165: VecPlaceArray(ctx->y2,py+m);
166: VecCopy(ctx->x1,ctx->y1);
167: MatMult(ctx->A,ctx->x2,ctx->w);
168: MatMult(ctx->AT,ctx->w,ctx->y2);
169: VecResetArray(ctx->x1);
170: VecResetArray(ctx->x2);
171: VecResetArray(ctx->y1);
172: VecResetArray(ctx->y2);
173: VecRestoreArrayRead(x,&px);
174: VecRestoreArrayWrite(y,&py);
175: return(0);
176: }
178: static PetscErrorCode MatGetDiagonal_ECross(Mat B,Vec d)
179: {
180: PetscErrorCode ierr;
181: SVD_CYCLIC_SHELL *ctx;
182: PetscScalar *pd;
183: PetscMPIInt len;
184: PetscInt mn,m,n,N,i,j,start,end,ncols;
185: PetscScalar *work1,*work2,*diag;
186: const PetscInt *cols;
187: const PetscScalar *vals;
190: MatShellGetContext(B,&ctx);
191: MatGetLocalSize(ctx->A,NULL,&n);
192: VecGetLocalSize(d,&mn);
193: m = mn-n;
194: VecGetArrayWrite(d,&pd);
195: VecPlaceArray(ctx->y1,pd);
196: VecSet(ctx->y1,1.0);
197: VecResetArray(ctx->y1);
198: VecPlaceArray(ctx->y2,pd+m);
199: if (!ctx->diag) {
200: /* compute diagonal from rows and store in ctx->diag */
201: VecDuplicate(ctx->y2,&ctx->diag);
202: MatGetSize(ctx->A,NULL,&N);
203: PetscCalloc2(N,&work1,N,&work2);
204: if (ctx->swapped) {
205: MatGetOwnershipRange(ctx->AT,&start,&end);
206: for (i=start;i<end;i++) {
207: MatGetRow(ctx->AT,i,&ncols,NULL,&vals);
208: for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
209: MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals);
210: }
211: } else {
212: MatGetOwnershipRange(ctx->A,&start,&end);
213: for (i=start;i<end;i++) {
214: MatGetRow(ctx->A,i,&ncols,&cols,&vals);
215: for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
216: MatRestoreRow(ctx->A,i,&ncols,&cols,&vals);
217: }
218: }
219: PetscMPIIntCast(N,&len);
220: MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B));
221: VecGetOwnershipRange(ctx->diag,&start,&end);
222: VecGetArrayWrite(ctx->diag,&diag);
223: for (i=start;i<end;i++) diag[i-start] = work2[i];
224: VecRestoreArrayWrite(ctx->diag,&diag);
225: PetscFree2(work1,work2);
226: }
227: VecCopy(ctx->diag,ctx->y2);
228: VecResetArray(ctx->y2);
229: VecRestoreArrayWrite(d,&pd);
230: return(0);
231: }
233: static PetscErrorCode MatDestroy_ECross(Mat B)
234: {
235: PetscErrorCode ierr;
236: SVD_CYCLIC_SHELL *ctx;
239: MatShellGetContext(B,&ctx);
240: VecDestroy(&ctx->x1);
241: VecDestroy(&ctx->x2);
242: VecDestroy(&ctx->y1);
243: VecDestroy(&ctx->y2);
244: VecDestroy(&ctx->diag);
245: VecDestroy(&ctx->w);
246: PetscFree(ctx);
247: return(0);
248: }
250: /*
251: Builds extended cross product matrix C = | I_m 0 |
252: | 0 AT*A |
253: t is an auxiliary Vec used to take the dimensions of the upper block
254: */
255: static PetscErrorCode SVDCyclicGetECrossMat(SVD svd,Mat A,Mat AT,Mat *C,Vec t)
256: {
257: PetscErrorCode ierr;
258: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
259: SVD_CYCLIC_SHELL *ctx;
260: PetscInt i,M,N,m,n,Istart,Iend;
261: VecType vtype;
262: Mat Id,Zm,Zn,ATA;
263: #if defined(PETSC_HAVE_CUDA)
264: PetscBool cuda;
265: #endif
268: MatGetSize(A,NULL,&N);
269: MatGetLocalSize(A,NULL,&n);
270: VecGetSize(t,&M);
271: VecGetLocalSize(t,&m);
273: if (cyclic->explicitmatrix) {
274: if (!svd->expltrans) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
275: MatCreate(PetscObjectComm((PetscObject)svd),&Id);
276: MatSetSizes(Id,m,m,M,M);
277: MatSetFromOptions(Id);
278: MatSetUp(Id);
279: MatGetOwnershipRange(Id,&Istart,&Iend);
280: for (i=Istart;i<Iend;i++) {
281: MatSetValue(Id,i,i,1.0,INSERT_VALUES);
282: }
283: MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
284: MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
285: MatCreate(PetscObjectComm((PetscObject)svd),&Zm);
286: MatSetSizes(Zm,m,n,M,N);
287: MatSetFromOptions(Zm);
288: MatSetUp(Zm);
289: MatGetOwnershipRange(Zm,&Istart,&Iend);
290: for (i=Istart;i<Iend;i++) {
291: if (i<N) { MatSetValue(Zm,i,i,0.0,INSERT_VALUES); }
292: }
293: MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY);
294: MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY);
295: MatCreate(PetscObjectComm((PetscObject)svd),&Zn);
296: MatSetSizes(Zn,n,m,N,M);
297: MatSetFromOptions(Zn);
298: MatSetUp(Zn);
299: MatGetOwnershipRange(Zn,&Istart,&Iend);
300: for (i=Istart;i<Iend;i++) {
301: if (i<m) { MatSetValue(Zn,i,i,0.0,INSERT_VALUES); }
302: }
303: MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY);
304: MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY);
305: MatProductCreate(AT,A,NULL,&ATA);
306: MatProductSetType(ATA,MATPRODUCT_AB);
307: MatProductSetFromOptions(ATA);
308: MatProductSymbolic(ATA);
309: MatProductNumeric(ATA);
310: MatCreateTile(1.0,Id,1.0,Zm,1.0,Zn,1.0,ATA,C);
311: MatDestroy(&Id);
312: MatDestroy(&Zm);
313: MatDestroy(&Zn);
314: MatDestroy(&ATA);
315: } else {
316: PetscNew(&ctx);
317: ctx->A = A;
318: ctx->AT = AT;
319: ctx->swapped = svd->swapped;
320: VecDuplicateEmpty(t,&ctx->x1);
321: VecDuplicateEmpty(t,&ctx->y1);
322: MatCreateVecsEmpty(A,&ctx->x2,NULL);
323: MatCreateVecsEmpty(A,&ctx->y2,NULL);
324: MatCreateVecs(A,NULL,&ctx->w);
325: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x1);
326: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->x2);
327: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y1);
328: PetscLogObjectParent((PetscObject)svd,(PetscObject)ctx->y2);
329: MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C);
330: MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ECross);
331: MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_ECross);
332: #if defined(PETSC_HAVE_CUDA)
333: PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&cuda,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,"");
334: if (cuda) {
335: MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross_CUDA);
336: } else
337: #endif
338: {
339: MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross);
340: }
341: MatGetVecType(A,&vtype);
342: MatSetVecType(*C,vtype);
343: }
344: PetscLogObjectParent((PetscObject)svd,(PetscObject)*C);
345: return(0);
346: }
348: PetscErrorCode SVDSetUp_Cyclic(SVD svd)
349: {
350: PetscErrorCode ierr;
351: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
352: PetscInt M,N,m,n,i,isl;
353: const PetscScalar *isa;
354: PetscScalar *va;
355: PetscBool trackall,issinv;
356: Vec v,t;
357: ST st;
360: MatGetSize(svd->A,&M,&N);
361: MatGetLocalSize(svd->A,&m,&n);
362: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
363: MatDestroy(&cyclic->C);
364: MatDestroy(&cyclic->D);
365: if (svd->isgeneralized) {
366: if (svd->which==SVD_SMALLEST) { /* alternative pencil */
367: MatCreateVecs(svd->B,NULL,&t);
368: SVDCyclicGetCyclicMat(svd,svd->B,svd->BT,&cyclic->C);
369: SVDCyclicGetECrossMat(svd,svd->A,svd->AT,&cyclic->D,t);
370: } else {
371: MatCreateVecs(svd->A,NULL,&t);
372: SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C);
373: SVDCyclicGetECrossMat(svd,svd->B,svd->BT,&cyclic->D,t);
374: }
375: VecDestroy(&t);
376: EPSSetOperators(cyclic->eps,cyclic->C,cyclic->D);
377: EPSSetProblemType(cyclic->eps,EPS_GHEP);
378: } else {
379: SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C);
380: EPSSetOperators(cyclic->eps,cyclic->C,NULL);
381: EPSSetProblemType(cyclic->eps,EPS_HEP);
382: }
383: if (!cyclic->usereps) {
384: if (svd->which == SVD_LARGEST) {
385: EPSGetST(cyclic->eps,&st);
386: PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv);
387: if (issinv) {
388: EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE);
389: } else {
390: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
391: }
392: } else {
393: if (svd->isgeneralized) { /* computes sigma^{-1} via alternative pencil */
394: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
395: } else {
396: EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL);
397: EPSSetTarget(cyclic->eps,0.0);
398: }
399: }
400: EPSSetDimensions(cyclic->eps,svd->nsv,svd->ncv,svd->mpd);
401: EPSSetTolerances(cyclic->eps,svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it);
402: switch (svd->conv) {
403: case SVD_CONV_ABS:
404: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_ABS);break;
405: case SVD_CONV_REL:
406: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_REL);break;
407: case SVD_CONV_NORM:
408: EPSSetConvergenceTest(cyclic->eps,EPS_CONV_NORM);break;
409: case SVD_CONV_MAXIT:
410: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
411: case SVD_CONV_USER:
412: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
413: }
414: }
415: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
416: /* Transfer the trackall option from svd to eps */
417: SVDGetTrackAll(svd,&trackall);
418: EPSSetTrackAll(cyclic->eps,trackall);
419: /* Transfer the initial subspace from svd to eps */
420: if (svd->nini<0 || svd->ninil<0) {
421: for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
422: MatCreateVecs(cyclic->C,&v,NULL);
423: VecGetArrayWrite(v,&va);
424: if (i<-svd->ninil) {
425: VecGetSize(svd->ISL[i],&isl);
426: if (isl!=m) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
427: VecGetArrayRead(svd->ISL[i],&isa);
428: PetscArraycpy(va,isa,m);
429: VecRestoreArrayRead(svd->IS[i],&isa);
430: } else {
431: PetscArrayzero(&va,m);
432: }
433: if (i<-svd->nini) {
434: VecGetSize(svd->IS[i],&isl);
435: if (isl!=n) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
436: VecGetArrayRead(svd->IS[i],&isa);
437: PetscArraycpy(va+m,isa,n);
438: VecRestoreArrayRead(svd->IS[i],&isa);
439: } else {
440: PetscArrayzero(va+m,n);
441: }
442: VecRestoreArrayWrite(v,&va);
443: VecDestroy(&svd->IS[i]);
444: svd->IS[i] = v;
445: }
446: svd->nini = PetscMin(svd->nini,svd->ninil);
447: EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS);
448: SlepcBasisDestroy_Private(&svd->nini,&svd->IS);
449: SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL);
450: }
451: EPSSetUp(cyclic->eps);
452: EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd);
453: svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
454: EPSGetTolerances(cyclic->eps,NULL,&svd->max_it);
455: if (svd->tol==PETSC_DEFAULT) svd->tol = SLEPC_DEFAULT_TOL;
457: svd->leftbasis = PETSC_TRUE;
458: SVDAllocateSolution(svd,0);
459: return(0);
460: }
462: PetscErrorCode SVDSolve_Cyclic(SVD svd)
463: {
465: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
466: PetscInt i,j,nconv;
467: PetscScalar sigma;
470: EPSSolve(cyclic->eps);
471: EPSGetConverged(cyclic->eps,&nconv);
472: EPSGetIterationNumber(cyclic->eps,&svd->its);
473: EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason);
474: for (i=0,j=0;i<nconv;i++) {
475: EPSGetEigenvalue(cyclic->eps,i,&sigma,NULL);
476: if (PetscRealPart(sigma) > 0.0) {
477: if (svd->isgeneralized && svd->which==SVD_SMALLEST) svd->sigma[j] = 1.0/PetscRealPart(sigma);
478: else svd->sigma[j] = PetscRealPart(sigma);
479: j++;
480: }
481: }
482: svd->nconv = j;
483: return(0);
484: }
486: PetscErrorCode SVDComputeVectors_Cyclic(SVD svd)
487: {
488: PetscErrorCode ierr;
489: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
490: PetscInt i,j,m,p,nconv;
491: PetscScalar *dst,sigma;
492: const PetscScalar *src,*px;
493: Vec u,v,x,x1,x2,uv;
496: EPSGetConverged(cyclic->eps,&nconv);
497: MatCreateVecs(cyclic->C,&x,NULL);
498: MatGetLocalSize(svd->A,&m,NULL);
499: if (svd->isgeneralized && svd->which==SVD_SMALLEST) {
500: MatCreateVecsEmpty(svd->B,&x1,&x2);
501: } else {
502: MatCreateVecsEmpty(svd->A,&x2,&x1);
503: }
504: if (svd->isgeneralized) {
505: MatCreateVecs(svd->A,NULL,&u);
506: MatCreateVecs(svd->B,NULL,&v);
507: MatGetLocalSize(svd->B,&p,NULL);
508: }
509: for (i=0,j=0;i<nconv;i++) {
510: EPSGetEigenpair(cyclic->eps,i,&sigma,NULL,x,NULL);
511: if (PetscRealPart(sigma) > 0.0) {
512: if (svd->isgeneralized) {
513: if (svd->which==SVD_SMALLEST) {
514: /* evec_i = 1/sqrt(2)*[ v_i; w_i ], w_i = x_i/c_i */
515: VecGetArrayRead(x,&px);
516: VecPlaceArray(x2,px);
517: VecPlaceArray(x1,px+p);
518: VecCopy(x2,v);
519: VecScale(v,PETSC_SQRT2); /* v_i = sqrt(2)*evec_i_1 */
520: VecScale(x1,PETSC_SQRT2); /* w_i = sqrt(2)*evec_i_2 */
521: MatMult(svd->A,x1,u); /* A*w_i = u_i */
522: VecScale(x1,1.0/PetscSqrtScalar(1.0+sigma*sigma)); /* x_i = w_i*c_i */
523: BVInsertVec(svd->V,j,x1);
524: VecResetArray(x2);
525: VecResetArray(x1);
526: VecRestoreArrayRead(x,&px);
527: } else {
528: /* evec_i = 1/sqrt(2)*[ u_i; w_i ], w_i = x_i/s_i */
529: VecGetArrayRead(x,&px);
530: VecPlaceArray(x1,px);
531: VecPlaceArray(x2,px+m);
532: VecCopy(x1,u);
533: VecScale(u,PETSC_SQRT2); /* u_i = sqrt(2)*evec_i_1 */
534: VecScale(x2,PETSC_SQRT2); /* w_i = sqrt(2)*evec_i_2 */
535: MatMult(svd->B,x2,v); /* B*w_i = v_i */
536: VecScale(x2,1.0/PetscSqrtScalar(1.0+sigma*sigma)); /* x_i = w_i*s_i */
537: BVInsertVec(svd->V,j,x2);
538: VecResetArray(x1);
539: VecResetArray(x2);
540: VecRestoreArrayRead(x,&px);
541: }
542: /* copy [u;v] to U[j] */
543: BVGetColumn(svd->U,j,&uv);
544: VecGetArrayWrite(uv,&dst);
545: VecGetArrayRead(u,&src);
546: PetscArraycpy(dst,src,m);
547: VecRestoreArrayRead(u,&src);
548: VecGetArrayRead(v,&src);
549: PetscArraycpy(dst+m,src,p);
550: VecRestoreArrayRead(v,&src);
551: VecRestoreArrayWrite(uv,&dst);
552: BVRestoreColumn(svd->U,j,&uv);
553: } else {
554: VecGetArrayRead(x,&px);
555: VecPlaceArray(x1,px);
556: VecPlaceArray(x2,px+m);
557: BVInsertVec(svd->U,j,x1);
558: BVScaleColumn(svd->U,j,PETSC_SQRT2);
559: BVInsertVec(svd->V,j,x2);
560: BVScaleColumn(svd->V,j,PETSC_SQRT2);
561: VecResetArray(x1);
562: VecResetArray(x2);
563: VecRestoreArrayRead(x,&px);
564: }
565: j++;
566: }
567: }
568: VecDestroy(&x);
569: VecDestroy(&x1);
570: VecDestroy(&x2);
571: if (svd->isgeneralized) {
572: VecDestroy(&u);
573: VecDestroy(&v);
574: }
575: return(0);
576: }
578: static PetscErrorCode EPSMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
579: {
580: PetscInt i,j;
581: SVD svd = (SVD)ctx;
582: PetscScalar er,ei;
586: nconv = 0;
587: for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
588: er = eigr[i]; ei = eigi[i];
589: STBackTransform(eps->st,1,&er,&ei);
590: if (PetscRealPart(er) > 0.0) {
591: svd->sigma[j] = PetscRealPart(er);
592: svd->errest[j] = errest[i];
593: if (errest[i] && errest[i] < svd->tol) nconv++;
594: j++;
595: }
596: }
597: nest = j;
598: SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest);
599: return(0);
600: }
602: PetscErrorCode SVDSetFromOptions_Cyclic(PetscOptionItems *PetscOptionsObject,SVD svd)
603: {
605: PetscBool set,val;
606: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
607: ST st;
610: PetscOptionsHead(PetscOptionsObject,"SVD Cyclic Options");
612: PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set);
613: if (set) { SVDCyclicSetExplicitMatrix(svd,val); }
615: PetscOptionsTail();
617: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
618: if (!cyclic->explicitmatrix && !cyclic->usereps) {
619: /* use as default an ST with shell matrix and Jacobi */
620: EPSGetST(cyclic->eps,&st);
621: STSetMatMode(st,ST_MATMODE_SHELL);
622: }
623: EPSSetFromOptions(cyclic->eps);
624: return(0);
625: }
627: static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmatrix)
628: {
629: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
632: if (cyclic->explicitmatrix != explicitmatrix) {
633: cyclic->explicitmatrix = explicitmatrix;
634: svd->state = SVD_STATE_INITIAL;
635: }
636: return(0);
637: }
639: /*@
640: SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
641: H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
643: Logically Collective on svd
645: Input Parameters:
646: + svd - singular value solver
647: - explicit - boolean flag indicating if H(A) is built explicitly
649: Options Database Key:
650: . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
652: Level: advanced
654: .seealso: SVDCyclicGetExplicitMatrix()
655: @*/
656: PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmatrix)
657: {
663: PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmatrix));
664: return(0);
665: }
667: static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmatrix)
668: {
669: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
672: *explicitmatrix = cyclic->explicitmatrix;
673: return(0);
674: }
676: /*@
677: SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly.
679: Not Collective
681: Input Parameter:
682: . svd - singular value solver
684: Output Parameter:
685: . explicit - the mode flag
687: Level: advanced
689: .seealso: SVDCyclicSetExplicitMatrix()
690: @*/
691: PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmatrix)
692: {
698: PetscUseMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmatrix));
699: return(0);
700: }
702: static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
703: {
704: PetscErrorCode ierr;
705: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
708: PetscObjectReference((PetscObject)eps);
709: EPSDestroy(&cyclic->eps);
710: cyclic->eps = eps;
711: cyclic->usereps = PETSC_TRUE;
712: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
713: svd->state = SVD_STATE_INITIAL;
714: return(0);
715: }
717: /*@
718: SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
719: singular value solver.
721: Collective on svd
723: Input Parameters:
724: + svd - singular value solver
725: - eps - the eigensolver object
727: Level: advanced
729: .seealso: SVDCyclicGetEPS()
730: @*/
731: PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
732: {
739: PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
740: return(0);
741: }
743: static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
744: {
746: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
749: if (!cyclic->eps) {
750: EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps);
751: PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1);
752: EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix);
753: EPSAppendOptionsPrefix(cyclic->eps,"svd_cyclic_");
754: PetscLogObjectParent((PetscObject)svd,(PetscObject)cyclic->eps);
755: PetscObjectSetOptions((PetscObject)cyclic->eps,((PetscObject)svd)->options);
756: EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL);
757: EPSMonitorSet(cyclic->eps,EPSMonitor_Cyclic,svd,NULL);
758: }
759: *eps = cyclic->eps;
760: return(0);
761: }
763: /*@
764: SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
765: to the singular value solver.
767: Not Collective
769: Input Parameter:
770: . svd - singular value solver
772: Output Parameter:
773: . eps - the eigensolver object
775: Level: advanced
777: .seealso: SVDCyclicSetEPS()
778: @*/
779: PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
780: {
786: PetscUseMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
787: return(0);
788: }
790: PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
791: {
793: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
794: PetscBool isascii;
797: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
798: if (isascii) {
799: if (!cyclic->eps) { SVDCyclicGetEPS(svd,&cyclic->eps); }
800: PetscViewerASCIIPrintf(viewer," %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit");
801: PetscViewerASCIIPushTab(viewer);
802: EPSView(cyclic->eps,viewer);
803: PetscViewerASCIIPopTab(viewer);
804: }
805: return(0);
806: }
808: PetscErrorCode SVDReset_Cyclic(SVD svd)
809: {
811: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
814: EPSReset(cyclic->eps);
815: MatDestroy(&cyclic->C);
816: MatDestroy(&cyclic->D);
817: return(0);
818: }
820: PetscErrorCode SVDDestroy_Cyclic(SVD svd)
821: {
823: SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
826: EPSDestroy(&cyclic->eps);
827: PetscFree(svd->data);
828: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL);
829: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL);
830: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL);
831: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL);
832: return(0);
833: }
835: SLEPC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
836: {
838: SVD_CYCLIC *cyclic;
841: PetscNewLog(svd,&cyclic);
842: svd->data = (void*)cyclic;
843: svd->ops->solve = SVDSolve_Cyclic;
844: svd->ops->solveg = SVDSolve_Cyclic;
845: svd->ops->setup = SVDSetUp_Cyclic;
846: svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
847: svd->ops->destroy = SVDDestroy_Cyclic;
848: svd->ops->reset = SVDReset_Cyclic;
849: svd->ops->view = SVDView_Cyclic;
850: svd->ops->computevectors = SVDComputeVectors_Cyclic;
851: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic);
852: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic);
853: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic);
854: PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic);
855: return(0);
856: }