Actual source code: linear.c
slepc-3.6.1 2015-09-03
1: /*
2: Explicit linearization for polynomial eigenproblems.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
25: #include linearp.h
29: static PetscErrorCode MatMult_Linear_Shift(Mat M,Vec x,Vec y)
30: {
31: PetscErrorCode ierr;
32: PEP_LINEAR *ctx;
33: PEP pep;
34: const PetscScalar *px;
35: PetscScalar *py,a,sigma=0.0;
36: PetscInt nmat,deg,i,m;
37: Vec x1,x2,x3,y1,aux;
38: PetscReal *ca,*cb,*cg;
39: PetscBool flg;
42: MatShellGetContext(M,(void**)&ctx);
43: pep = ctx->pep;
44: STGetTransform(pep->st,&flg);
45: if (!flg) {
46: STGetShift(pep->st,&sigma);
47: }
48: nmat = pep->nmat;
49: deg = nmat-1;
50: m = pep->nloc;
51: ca = pep->pbc;
52: cb = pep->pbc+nmat;
53: cg = pep->pbc+2*nmat;
54: x1=ctx->w[0];x2=ctx->w[1];x3=ctx->w[2];y1=ctx->w[3];aux=ctx->w[4];
55:
56: VecSet(y,0.0);
57: VecGetArrayRead(x,&px);
58: VecGetArray(y,&py);
59: a = 1.0;
61: /* first block */
62: VecPlaceArray(x2,px);
63: VecPlaceArray(x3,px+m);
64: VecPlaceArray(y1,py);
65: VecAXPY(y1,cb[0]-sigma,x2);
66: VecAXPY(y1,ca[0],x3);
67: VecResetArray(x2);
68: VecResetArray(x3);
69: VecResetArray(y1);
71: /* inner blocks */
72: for (i=1;i<deg-1;i++) {
73: VecPlaceArray(x1,px+(i-1)*m);
74: VecPlaceArray(x2,px+i*m);
75: VecPlaceArray(x3,px+(i+1)*m);
76: VecPlaceArray(y1,py+i*m);
77: VecAXPY(y1,cg[i],x1);
78: VecAXPY(y1,cb[i]-sigma,x2);
79: VecAXPY(y1,ca[i],x3);
80: VecResetArray(x1);
81: VecResetArray(x2);
82: VecResetArray(x3);
83: VecResetArray(y1);
84: }
86: /* last block */
87: VecPlaceArray(y1,py+(deg-1)*m);
88: for (i=0;i<deg;i++) {
89: VecPlaceArray(x1,px+i*m);
90: STMatMult(pep->st,i,x1,aux);
91: VecAXPY(y1,a,aux);
92: VecResetArray(x1);
93: a *= pep->sfactor;
94: }
95: VecCopy(y1,aux);
96: STMatSolve(pep->st,aux,y1);
97: VecScale(y1,-ca[deg-1]/a);
98: VecPlaceArray(x1,px+(deg-2)*m);
99: VecPlaceArray(x2,px+(deg-1)*m);
100: VecAXPY(y1,cg[deg-1],x1);
101: VecAXPY(y1,cb[deg-1]-sigma,x2);
102: VecResetArray(x1);
103: VecResetArray(x2);
104: VecResetArray(y1);
106: VecRestoreArrayRead(x,&px);
107: VecRestoreArray(y,&py);
108: return(0);
109: }
113: static PetscErrorCode MatMult_Linear_Sinvert(Mat M,Vec x,Vec y)
114: {
115: PetscErrorCode ierr;
116: PEP_LINEAR *ctx;
117: PEP pep;
118: const PetscScalar *px;
119: PetscScalar *py,a,sigma,t=1.0,tp=0.0,tt;
120: PetscInt nmat,deg,i,m;
121: Vec x1,y1,y2,y3,aux,aux2;
122: PetscReal *ca,*cb,*cg;
125: MatShellGetContext(M,(void**)&ctx);
126: pep = ctx->pep;
127: nmat = pep->nmat;
128: deg = nmat-1;
129: m = pep->nloc;
130: ca = pep->pbc;
131: cb = pep->pbc+nmat;
132: cg = pep->pbc+2*nmat;
133: x1=ctx->w[0];y1=ctx->w[1];y2=ctx->w[2];y3=ctx->w[3];aux=ctx->w[4];aux2=ctx->w[5];
134: EPSGetTarget(ctx->eps,&sigma);
135: VecSet(y,0.0);
136: VecGetArrayRead(x,&px);
137: VecGetArray(y,&py);
138: a = pep->sfactor;
140: /* first block */
141: VecPlaceArray(x1,px);
142: VecPlaceArray(y1,py+m);
143: VecCopy(x1,y1);
144: VecScale(y1,1.0/ca[0]);
145: VecResetArray(x1);
146: VecResetArray(y1);
148: /* second block */
149: if (deg>2) {
150: VecPlaceArray(x1,px+m);
151: VecPlaceArray(y1,py+m);
152: VecPlaceArray(y2,py+2*m);
153: VecCopy(x1,y2);
154: VecAXPY(y2,sigma-cb[1],y1);
155: VecScale(y2,1.0/ca[1]);
156: VecResetArray(x1);
157: VecResetArray(y1);
158: VecResetArray(y2);
159: }
161: /* inner blocks */
162: for (i=2;i<deg-1;i++) {
163: VecPlaceArray(x1,px+i*m);
164: VecPlaceArray(y1,py+(i-1)*m);
165: VecPlaceArray(y2,py+i*m);
166: VecPlaceArray(y3,py+(i+1)*m);
167: VecCopy(x1,y3);
168: VecAXPY(y3,sigma-cb[i],y2);
169: VecAXPY(y3,-cg[i],y1);
170: VecScale(y3,1.0/ca[i]);
171: VecResetArray(x1);
172: VecResetArray(y1);
173: VecResetArray(y2);
174: VecResetArray(y3);
175: }
177: /* last block */
178: VecPlaceArray(y1,py);
179: for (i=0;i<deg-2;i++) {
180: VecPlaceArray(y2,py+(i+1)*m);
181: STMatMult(pep->st,i+1,y2,aux);
182: VecAXPY(y1,a,aux);
183: VecResetArray(y2);
184: a *= pep->sfactor;
185: }
186: i = deg-2;
187: VecPlaceArray(y2,py+(i+1)*m);
188: VecPlaceArray(y3,py+i*m);
189: VecCopy(y2,aux2);
190: VecAXPY(aux2,cg[i+1]/ca[i+1],y3);
191: STMatMult(pep->st,i+1,aux2,aux);
192: VecAXPY(y1,a,aux);
193: VecResetArray(y2);
194: VecResetArray(y3);
195: a *= pep->sfactor;
196: i = deg-1;
197: VecPlaceArray(x1,px+i*m);
198: VecPlaceArray(y3,py+i*m);
199: VecCopy(x1,aux2);
200: VecAXPY(aux2,sigma-cb[i],y3);
201: VecScale(aux2,1.0/ca[i]);
202: STMatMult(pep->st,i+1,aux2,aux);
203: VecAXPY(y1,a,aux);
204: VecResetArray(x1);
205: VecResetArray(y3);
207: VecCopy(y1,aux);
208: STMatSolve(pep->st,aux,y1);
209: VecScale(y1,-1.0);
211: /* final update */
212: for (i=1;i<deg;i++) {
213: VecPlaceArray(y2,py+i*m);
214: tt = t;
215: t = ((sigma-cb[i-1])*t-cg[i-1]*tp)/ca[i-1]; /* i-th basis polynomial */
216: tp = tt;
217: VecAXPY(y2,t,y1);
218: VecResetArray(y2);
219: }
220: VecResetArray(y1);
222: VecRestoreArrayRead(x,&px);
223: VecRestoreArray(y,&py);
224: return(0);
225: }
229: static PetscErrorCode BackTransform_Linear(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
230: {
232: PEP_LINEAR *ctx;
233: ST stctx;
236: STShellGetContext(st,(void**)&ctx);
237: PEPGetST(ctx->pep,&stctx);
238: STBackTransform(stctx,n,eigr,eigi);
239: return(0);
240: }
244: static PetscErrorCode Apply_Linear(ST st,Vec x,Vec y)
245: {
247: PEP_LINEAR *ctx;
250: STShellGetContext(st,(void**)&ctx);
251: MatMult(ctx->A,x,y);
252: return(0);
253: }
257: PetscErrorCode PEPSetUp_Linear(PEP pep)
258: {
260: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
261: ST st;
262: PetscInt i=0,deg=pep->nmat-1;
263: EPSWhich which;
264: EPSProblemType ptype;
265: PetscBool trackall,istrivial,transf,sinv,ks;
266: PetscScalar sigma,*epsarray,*peparray;
267: Vec veps;
268: /* function tables */
269: PetscErrorCode (*fcreate[][2])(MPI_Comm,PEP_LINEAR*,Mat*) = {
270: { MatCreateExplicit_Linear_N1A, MatCreateExplicit_Linear_N1B }, /* N1 */
271: { MatCreateExplicit_Linear_N2A, MatCreateExplicit_Linear_N2B }, /* N2 */
272: { MatCreateExplicit_Linear_S1A, MatCreateExplicit_Linear_S1B }, /* S1 */
273: { MatCreateExplicit_Linear_S2A, MatCreateExplicit_Linear_S2B }, /* S2 */
274: { MatCreateExplicit_Linear_H1A, MatCreateExplicit_Linear_H1B }, /* H1 */
275: { MatCreateExplicit_Linear_H2A, MatCreateExplicit_Linear_H2B } /* H2 */
276: };
279: pep->lineariz = PETSC_TRUE;
280: if (!ctx->cform) ctx->cform = 1;
281: STGetTransform(pep->st,&transf);
282: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
283: if (!pep->which) {
284: if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
285: else pep->which = PEP_LARGEST_MAGNITUDE;
286: }
287: STSetUp(pep->st);
288: if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
289: EPSGetST(ctx->eps,&st);
290: if (!transf) { EPSSetTarget(ctx->eps,pep->target); }
291: if (sinv && !transf) { STSetDefaultShift(st,pep->target); }
292: /* compute scale factor if not set by user */
293: PEPComputeScaleFactor(pep);
295: if (ctx->explicitmatrix) {
296: if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option is not implemented with st-tranform flag active");
297: if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option only available for quadratic problems");
298: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option not implemented for non-monomial bases");
299: if (pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Diagonal scaling not allowed in PEPLINEAR with explicit matrices");
300: if (sinv && !transf) { STSetType(st,STSINVERT); }
301: RGSetScale(pep->rg,pep->sfactor);
302: STGetTOperators(pep->st,0,&ctx->K);
303: STGetTOperators(pep->st,1,&ctx->C);
304: STGetTOperators(pep->st,2,&ctx->M);
305: ctx->sfactor = pep->sfactor;
306: ctx->dsfactor = pep->dsfactor;
307:
308: MatDestroy(&ctx->A);
309: MatDestroy(&ctx->B);
310: VecDestroy(&ctx->w[0]);
311: VecDestroy(&ctx->w[1]);
312: VecDestroy(&ctx->w[2]);
313: VecDestroy(&ctx->w[3]);
314:
315: switch (pep->problem_type) {
316: case PEP_GENERAL: i = 0; break;
317: case PEP_HERMITIAN: i = 2; break;
318: case PEP_GYROSCOPIC: i = 4; break;
319: default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->problem_type");
320: }
321: i += ctx->cform-1;
323: (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
324: (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);
325: PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
326: PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);
328: } else { /* implicit matrix */
329: if (pep->problem_type!=PEP_GENERAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Must use the explicit matrix option");
330: if (!((PetscObject)(ctx->eps))->type_name) {
331: EPSSetType(ctx->eps,EPSKRYLOVSCHUR);
332: } else {
333: PetscObjectTypeCompare((PetscObject)ctx->eps,EPSKRYLOVSCHUR,&ks);
334: if (!ks) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option only implemented for Krylov-Schur");
335: }
336: if (ctx->cform!=1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option not available for 2nd companion form");
337: STSetType(st,STSHELL);
338: STShellSetContext(st,(PetscObject)ctx);
339: if (!transf) { STShellSetBackTransform(st,BackTransform_Linear); }
340: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[0]);
341: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[1]);
342: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[2]);
343: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[3]);
344: MatCreateVecs(pep->A[0],&ctx->w[4],NULL);
345: MatCreateVecs(pep->A[0],&ctx->w[5],NULL);
346: PetscLogObjectParents(pep,6,ctx->w);
347: MatCreateShell(PetscObjectComm((PetscObject)pep),deg*pep->nloc,deg*pep->nloc,deg*pep->n,deg*pep->n,ctx,&ctx->A);
348: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
349: if (sinv && !transf) {
350: MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Sinvert);
351: } else {
352: MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Shift);
353: }
354: STShellSetApply(st,Apply_Linear);
355: PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
356: ctx->pep = pep;
358: PEPBasisCoefficients(pep,pep->pbc);
359: if (!transf) {
360: PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
361: if (sinv) {
362: PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
363: } else {
364: for (i=0;i<deg;i++) pep->solvematcoeffs[i] = 0.0;
365: pep->solvematcoeffs[deg] = 1.0;
366: }
367: STScaleShift(pep->st,1.0/pep->sfactor);
368: RGSetScale(pep->rg,pep->sfactor);
369: }
370: if (pep->sfactor!=1.0) {
371: for (i=0;i<pep->nmat;i++) {
372: pep->pbc[pep->nmat+i] /= pep->sfactor;
373: pep->pbc[2*pep->nmat+i] /= pep->sfactor*pep->sfactor;
374: }
375: }
376: }
378: EPSSetOperators(ctx->eps,ctx->A,ctx->B);
379: EPSGetProblemType(ctx->eps,&ptype);
380: if (!ptype) {
381: if (ctx->explicitmatrix) {
382: EPSSetProblemType(ctx->eps,EPS_GNHEP);
383: } else {
384: EPSSetProblemType(ctx->eps,EPS_NHEP);
385: }
386: }
387: if (transf) which = EPS_LARGEST_MAGNITUDE;
388: else {
389: switch (pep->which) {
390: case PEP_LARGEST_MAGNITUDE: which = EPS_LARGEST_MAGNITUDE; break;
391: case PEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
392: case PEP_LARGEST_REAL: which = EPS_LARGEST_REAL; break;
393: case PEP_SMALLEST_REAL: which = EPS_SMALLEST_REAL; break;
394: case PEP_LARGEST_IMAGINARY: which = EPS_LARGEST_IMAGINARY; break;
395: case PEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
396: case PEP_TARGET_MAGNITUDE: which = EPS_TARGET_MAGNITUDE; break;
397: case PEP_TARGET_REAL: which = EPS_TARGET_REAL; break;
398: case PEP_TARGET_IMAGINARY: which = EPS_TARGET_IMAGINARY; break;
399: case PEP_WHICH_USER: which = EPS_WHICH_USER;
400: EPSSetEigenvalueComparison(ctx->eps,pep->sc->comparison,pep->sc->comparisonctx);
401: break;
402: default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of which");
403: }
404: }
405: EPSSetWhichEigenpairs(ctx->eps,which);
407: EPSSetDimensions(ctx->eps,pep->nev,pep->ncv?pep->ncv:PETSC_DEFAULT,pep->mpd?pep->mpd:PETSC_DEFAULT);
408: EPSSetTolerances(ctx->eps,pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,pep->max_it?pep->max_it:PETSC_DEFAULT);
409: RGIsTrivial(pep->rg,&istrivial);
410: if (!istrivial) {
411: if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEPLINEAR does not support a nontrivial region with st-transform");
412: EPSSetRG(ctx->eps,pep->rg);
413: }
414: /* Transfer the trackall option from pep to eps */
415: PEPGetTrackAll(pep,&trackall);
416: EPSSetTrackAll(ctx->eps,trackall);
418: /* temporary change of target */
419: if (pep->sfactor!=1.0) {
420: EPSGetTarget(ctx->eps,&sigma);
421: EPSSetTarget(ctx->eps,sigma/pep->sfactor);
422: }
424: /* process initial vector */
425: if (pep->nini<=-deg) {
426: VecCreateMPI(PetscObjectComm((PetscObject)ctx->eps),deg*pep->nloc,deg*pep->n,&veps);
427: VecGetArray(veps,&epsarray);
428: for (i=0;i<deg;i++) {
429: VecGetArray(pep->IS[i],&peparray);
430: PetscMemcpy(epsarray+i*pep->nloc,peparray,pep->nloc*sizeof(PetscScalar));
431: VecRestoreArray(pep->IS[i],&peparray);
432: }
433: VecRestoreArray(veps,&epsarray);
434: EPSSetInitialSpace(ctx->eps,1,&veps);
435: VecDestroy(&veps);
436: }
437: if (pep->nini<0) {
438: SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
439: }
441: EPSSetUp(ctx->eps);
442: EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
443: EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
444: if (pep->nini>0) { PetscInfo(pep,"Ignoring initial vectors\n"); }
445: PEPAllocateSolution(pep,0);
446: return(0);
447: }
451: /*
452: PEPLinearExtract_Residual - Auxiliary routine that copies the solution of the
453: linear eigenproblem to the PEP object. The eigenvector of the generalized
454: problem is supposed to be
455: z = [ x ]
456: [ l*x ]
457: The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
458: computed residual norm.
459: Finally, x is normalized so that ||x||_2 = 1.
460: */
461: static PetscErrorCode PEPLinearExtract_Residual(PEP pep,EPS eps)
462: {
463: PetscErrorCode ierr;
464: PetscInt i,k;
465: const PetscScalar *px;
466: PetscScalar *er=pep->eigr,*ei=pep->eigi;
467: PetscReal rn1,rn2;
468: Vec xr,xi=NULL,wr;
469: Mat A;
470: #if !defined(PETSC_USE_COMPLEX)
471: Vec wi;
472: const PetscScalar *py;
473: #endif
476: #if defined(PETSC_USE_COMPLEX)
477: PEPSetWorkVecs(pep,2);
478: #else
479: PEPSetWorkVecs(pep,4);
480: #endif
481: EPSGetOperators(eps,&A,NULL);
482: MatCreateVecs(A,&xr,NULL);
483: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&wr);
484: #if !defined(PETSC_USE_COMPLEX)
485: VecDuplicate(xr,&xi);
486: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&wi);
487: #endif
488: for (i=0;i<pep->nconv;i++) {
489: EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
490: #if !defined(PETSC_USE_COMPLEX)
491: if (ei[i]!=0.0) { /* complex conjugate pair */
492: VecGetArrayRead(xr,&px);
493: VecGetArrayRead(xi,&py);
494: VecPlaceArray(wr,px);
495: VecPlaceArray(wi,py);
496: SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
497: PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn1);
498: BVInsertVec(pep->V,i,wr);
499: BVInsertVec(pep->V,i+1,wi);
500: for (k=1;k<pep->nmat-1;k++) {
501: VecResetArray(wr);
502: VecResetArray(wi);
503: VecPlaceArray(wr,px+k*pep->nloc);
504: VecPlaceArray(wi,py+k*pep->nloc);
505: SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
506: PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn2);
507: if (rn1>rn2) {
508: BVInsertVec(pep->V,i,wr);
509: BVInsertVec(pep->V,i+1,wi);
510: rn1 = rn2;
511: }
512: }
513: VecResetArray(wr);
514: VecResetArray(wi);
515: VecRestoreArrayRead(xr,&px);
516: VecRestoreArrayRead(xi,&py);
517: i++;
518: } else /* real eigenvalue */
519: #endif
520: {
521: VecGetArrayRead(xr,&px);
522: VecPlaceArray(wr,px);
523: SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
524: PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn1);
525: BVInsertVec(pep->V,i,wr);
526: for (k=1;k<pep->nmat-1;k++) {
527: VecResetArray(wr);
528: VecPlaceArray(wr,px+k*pep->nloc);
529: SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
530: PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn2);
531: if (rn1>rn2) {
532: BVInsertVec(pep->V,i,wr);
533: rn1 = rn2;
534: }
535: }
536: VecResetArray(wr);
537: VecRestoreArrayRead(xr,&px);
538: }
539: }
540: VecDestroy(&wr);
541: VecDestroy(&xr);
542: #if !defined(PETSC_USE_COMPLEX)
543: VecDestroy(&wi);
544: VecDestroy(&xi);
545: #endif
546: return(0);
547: }
551: /*
552: PEPLinearExtract_None - Same as PEPLinearExtract_Norm but always takes
553: the first block.
554: */
555: static PetscErrorCode PEPLinearExtract_None(PEP pep,EPS eps)
556: {
557: PetscErrorCode ierr;
558: PetscInt i;
559: const PetscScalar *px;
560: Mat A;
561: Vec xr,xi,w;
562: #if !defined(PETSC_USE_COMPLEX)
563: PetscScalar *ei=pep->eigi;
564: #endif
567: EPSGetOperators(eps,&A,NULL);
568: MatCreateVecs(A,&xr,NULL);
569: VecDuplicate(xr,&xi);
570: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&w);
571: for (i=0;i<pep->nconv;i++) {
572: EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
573: #if !defined(PETSC_USE_COMPLEX)
574: if (ei[i]!=0.0) { /* complex conjugate pair */
575: VecGetArrayRead(xr,&px);
576: VecPlaceArray(w,px);
577: BVInsertVec(pep->V,i,w);
578: VecResetArray(w);
579: VecRestoreArrayRead(xr,&px);
580: VecGetArrayRead(xi,&px);
581: VecPlaceArray(w,px);
582: BVInsertVec(pep->V,i+1,w);
583: VecResetArray(w);
584: VecRestoreArrayRead(xi,&px);
585: i++;
586: } else /* real eigenvalue */
587: #endif
588: {
589: VecGetArrayRead(xr,&px);
590: VecPlaceArray(w,px);
591: BVInsertVec(pep->V,i,w);
592: VecResetArray(w);
593: VecRestoreArrayRead(xr,&px);
594: }
595: }
596: VecDestroy(&w);
597: VecDestroy(&xr);
598: VecDestroy(&xi);
599: return(0);
600: }
604: /*
605: PEPLinearExtract_Norm - Auxiliary routine that copies the solution of the
606: linear eigenproblem to the PEP object. The eigenvector of the generalized
607: problem is supposed to be
608: z = [ x ]
609: [ l*x ]
610: If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
611: Finally, x is normalized so that ||x||_2 = 1.
612: */
613: static PetscErrorCode PEPLinearExtract_Norm(PEP pep,EPS eps)
614: {
615: PetscErrorCode ierr;
616: PetscInt i,offset;
617: const PetscScalar *px;
618: PetscScalar *er=pep->eigr;
619: Mat A;
620: Vec xr,xi=NULL,w;
621: #if !defined(PETSC_USE_COMPLEX)
622: PetscScalar *ei=pep->eigi;
623: #endif
626: EPSGetOperators(eps,&A,NULL);
627: MatCreateVecs(A,&xr,NULL);
628: #if !defined(PETSC_USE_COMPLEX)
629: VecDuplicate(xr,&xi);
630: #endif
631: VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&w);
632: for (i=0;i<pep->nconv;i++) {
633: EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
634: if (SlepcAbsEigenvalue(er[i],ei[i])>1.0) offset = (pep->nmat-2)*pep->nloc;
635: else offset = 0;
636: #if !defined(PETSC_USE_COMPLEX)
637: if (ei[i]!=0.0) { /* complex conjugate pair */
638: VecGetArrayRead(xr,&px);
639: VecPlaceArray(w,px+offset);
640: BVInsertVec(pep->V,i,w);
641: VecResetArray(w);
642: VecRestoreArrayRead(xr,&px);
643: VecGetArrayRead(xi,&px);
644: VecPlaceArray(w,px+offset);
645: BVInsertVec(pep->V,i+1,w);
646: VecResetArray(w);
647: VecRestoreArrayRead(xi,&px);
648: i++;
649: } else /* real eigenvalue */
650: #endif
651: {
652: VecGetArrayRead(xr,&px);
653: VecPlaceArray(w,px+offset);
654: BVInsertVec(pep->V,i,w);
655: VecResetArray(w);
656: VecRestoreArrayRead(xr,&px);
657: }
658: }
659: VecDestroy(&w);
660: VecDestroy(&xr);
661: #if !defined(PETSC_USE_COMPLEX)
662: VecDestroy(&xi);
663: #endif
664: return(0);
665: }
669: PetscErrorCode PEPExtractVectors_Linear(PEP pep)
670: {
672: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
673:
675: switch (pep->extract) {
676: case PEP_EXTRACT_NONE:
677: PEPLinearExtract_None(pep,ctx->eps);
678: break;
679: case PEP_EXTRACT_NORM:
680: PEPLinearExtract_Norm(pep,ctx->eps);
681: break;
682: case PEP_EXTRACT_RESIDUAL:
683: PEPLinearExtract_Residual(pep,ctx->eps);
684: break;
685: default:
686: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
687: }
688: return(0);
689: }
693: PetscErrorCode PEPSolve_Linear(PEP pep)
694: {
696: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
697: PetscScalar sigma;
698: PetscBool flg;
699: PetscInt i;
702: EPSSolve(ctx->eps);
703: EPSGetConverged(ctx->eps,&pep->nconv);
704: EPSGetIterationNumber(ctx->eps,&pep->its);
705: EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);
707: /* recover eigenvalues */
708: for (i=0;i<pep->nconv;i++) {
709: EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
710: pep->eigr[i] *= pep->sfactor;
711: pep->eigi[i] *= pep->sfactor;
712: }
714: /* restore target */
715: EPSGetTarget(ctx->eps,&sigma);
716: EPSSetTarget(ctx->eps,sigma*pep->sfactor);
718: STGetTransform(pep->st,&flg);
719: if (flg && pep->ops->backtransform) {
720: (*pep->ops->backtransform)(pep);
721: }
722: if (pep->sfactor!=1.0) {
723: /* Restore original values */
724: for (i=0;i<pep->nmat;i++){
725: pep->pbc[pep->nmat+i] *= pep->sfactor;
726: pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
727: }
728: if (!flg && !ctx->explicitmatrix) {
729: STScaleShift(pep->st,pep->sfactor);
730: }
731: RGSetScale(pep->rg,1.0);
732: }
733: return(0);
734: }
738: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
739: {
740: PetscInt i;
741: PEP pep = (PEP)ctx;
742: ST st;
746: for (i=0;i<PetscMin(nest,pep->ncv);i++) {
747: pep->eigr[i] = eigr[i];
748: pep->eigi[i] = eigi[i];
749: pep->errest[i] = errest[i];
750: }
751: EPSGetST(eps,&st);
752: STBackTransform(st,nest,pep->eigr,pep->eigi);
753: PEPMonitor(pep,its,nconv,pep->eigr,pep->eigi,pep->errest,nest);
754: return(0);
755: }
759: PetscErrorCode PEPSetFromOptions_Linear(PetscOptions *PetscOptionsObject,PEP pep)
760: {
762: PetscBool set,val;
763: PetscInt i;
764: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
767: PetscOptionsHead(PetscOptionsObject,"PEP Linear Options");
768: PetscOptionsInt("-pep_linear_cform","Number of the companion form","PEPLinearSetCompanionForm",ctx->cform,&i,&set);
769: if (set) {
770: PEPLinearSetCompanionForm(pep,i);
771: }
772: PetscOptionsBool("-pep_linear_explicitmatrix","Use explicit matrix in linearization","PEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
773: if (set) {
774: PEPLinearSetExplicitMatrix(pep,val);
775: }
776: if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
777: EPSSetFromOptions(ctx->eps);
778: PetscOptionsTail();
779: return(0);
780: }
784: static PetscErrorCode PEPLinearSetCompanionForm_Linear(PEP pep,PetscInt cform)
785: {
786: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
789: if (!cform) return(0);
790: if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
791: else {
792: if (cform!=1 && cform!=2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
793: ctx->cform = cform;
794: }
795: return(0);
796: }
800: /*@
801: PEPLinearSetCompanionForm - Choose between the two companion forms available
802: for the linearization of a quadratic eigenproblem.
804: Logically Collective on PEP
806: Input Parameters:
807: + pep - polynomial eigenvalue solver
808: - cform - 1 or 2 (first or second companion form)
810: Options Database Key:
811: . -pep_linear_cform <int> - Choose the companion form
813: Level: advanced
815: .seealso: PEPLinearGetCompanionForm()
816: @*/
817: PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform)
818: {
824: PetscTryMethod(pep,"PEPLinearSetCompanionForm_C",(PEP,PetscInt),(pep,cform));
825: return(0);
826: }
830: static PetscErrorCode PEPLinearGetCompanionForm_Linear(PEP pep,PetscInt *cform)
831: {
832: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
835: *cform = ctx->cform;
836: return(0);
837: }
841: /*@
842: PEPLinearGetCompanionForm - Returns the number of the companion form that
843: will be used for the linearization of a quadratic eigenproblem.
845: Not Collective
847: Input Parameter:
848: . pep - polynomial eigenvalue solver
850: Output Parameter:
851: . cform - the companion form number (1 or 2)
853: Level: advanced
855: .seealso: PEPLinearSetCompanionForm()
856: @*/
857: PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform)
858: {
864: PetscTryMethod(pep,"PEPLinearGetCompanionForm_C",(PEP,PetscInt*),(pep,cform));
865: return(0);
866: }
870: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
871: {
872: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
875: ctx->explicitmatrix = explicitmatrix;
876: return(0);
877: }
881: /*@
882: PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
883: linearization of the problem must be built explicitly.
885: Logically Collective on PEP
887: Input Parameters:
888: + pep - polynomial eigenvalue solver
889: - explicit - boolean flag indicating if the matrices are built explicitly
891: Options Database Key:
892: . -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag
894: Level: advanced
896: .seealso: PEPLinearGetExplicitMatrix()
897: @*/
898: PetscErrorCode PEPLinearSetExplicitMatrix(PEP pep,PetscBool explicitmatrix)
899: {
905: PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmatrix));
906: return(0);
907: }
911: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmatrix)
912: {
913: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
916: *explicitmatrix = ctx->explicitmatrix;
917: return(0);
918: }
922: /*@
923: PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
924: A and B for the linearization are built explicitly.
926: Not Collective
928: Input Parameter:
929: . pep - polynomial eigenvalue solver
931: Output Parameter:
932: . explicitmatrix - the mode flag
934: Level: advanced
936: .seealso: PEPLinearSetExplicitMatrix()
937: @*/
938: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmatrix)
939: {
945: PetscTryMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmatrix));
946: return(0);
947: }
951: static PetscErrorCode PEPLinearSetEPS_Linear(PEP pep,EPS eps)
952: {
954: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
957: PetscObjectReference((PetscObject)eps);
958: EPSDestroy(&ctx->eps);
959: ctx->eps = eps;
960: PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
961: pep->state = PEP_STATE_INITIAL;
962: return(0);
963: }
967: /*@
968: PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
969: polynomial eigenvalue solver.
971: Collective on PEP
973: Input Parameters:
974: + pep - polynomial eigenvalue solver
975: - eps - the eigensolver object
977: Level: advanced
979: .seealso: PEPLinearGetEPS()
980: @*/
981: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
982: {
989: PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
990: return(0);
991: }
995: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
996: {
998: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
999: ST st;
1002: if (!ctx->eps) {
1003: EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
1004: EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
1005: EPSAppendOptionsPrefix(ctx->eps,"pep_");
1006: EPSGetST(ctx->eps,&st);
1007: STSetOptionsPrefix(st,((PetscObject)ctx->eps)->prefix);
1008: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
1009: PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
1010: EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
1011: }
1012: *eps = ctx->eps;
1013: return(0);
1014: }
1018: /*@
1019: PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
1020: to the polynomial eigenvalue solver.
1022: Not Collective
1024: Input Parameter:
1025: . pep - polynomial eigenvalue solver
1027: Output Parameter:
1028: . eps - the eigensolver object
1030: Level: advanced
1032: .seealso: PEPLinearSetEPS()
1033: @*/
1034: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
1035: {
1041: PetscTryMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
1042: return(0);
1043: }
1047: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
1048: {
1050: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
1053: if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
1054: PetscViewerASCIIPrintf(viewer," Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
1055: PetscViewerASCIIPrintf(viewer," Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
1056: PetscViewerASCIIPushTab(viewer);
1057: EPSView(ctx->eps,viewer);
1058: PetscViewerASCIIPopTab(viewer);
1059: return(0);
1060: }
1064: PetscErrorCode PEPReset_Linear(PEP pep)
1065: {
1067: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
1070: if (!ctx->eps) { EPSReset(ctx->eps); }
1071: MatDestroy(&ctx->A);
1072: MatDestroy(&ctx->B);
1073: VecDestroy(&ctx->w[0]);
1074: VecDestroy(&ctx->w[1]);
1075: VecDestroy(&ctx->w[2]);
1076: VecDestroy(&ctx->w[3]);
1077: VecDestroy(&ctx->w[4]);
1078: VecDestroy(&ctx->w[5]);
1079: return(0);
1080: }
1084: PetscErrorCode PEPDestroy_Linear(PEP pep)
1085: {
1087: PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;
1090: EPSDestroy(&ctx->eps);
1091: PetscFree(pep->data);
1092: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",NULL);
1093: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",NULL);
1094: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
1095: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
1096: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
1097: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
1098: return(0);
1099: }
1103: PETSC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
1104: {
1106: PEP_LINEAR *ctx;
1109: PetscNewLog(pep,&ctx);
1110: ctx->explicitmatrix = PETSC_FALSE;
1111: pep->data = (void*)ctx;
1113: pep->ops->solve = PEPSolve_Linear;
1114: pep->ops->setup = PEPSetUp_Linear;
1115: pep->ops->setfromoptions = PEPSetFromOptions_Linear;
1116: pep->ops->destroy = PEPDestroy_Linear;
1117: pep->ops->reset = PEPReset_Linear;
1118: pep->ops->view = PEPView_Linear;
1119: pep->ops->backtransform = PEPBackTransform_Default;
1120: pep->ops->computevectors = PEPComputeVectors_Default;
1121: pep->ops->extractvectors = PEPExtractVectors_Linear;
1122: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",PEPLinearSetCompanionForm_Linear);
1123: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",PEPLinearGetCompanionForm_Linear);
1124: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
1125: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
1126: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
1127: PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
1128: return(0);
1129: }