Actual source code: peprefine.c
slepc-3.7.1 2016-05-27
1: /*
2: Newton refinement for PEP, simple version.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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>
25: #include <slepcblaslapack.h>
27: #define NREF_MAXIT 10
29: typedef struct {
30: VecScatter *scatter_id,nst;
31: Mat *A;
32: Vec nv,vg,v,w;
33: } PEPSimpNRefctx;
35: typedef struct {
36: Mat M1;
37: Vec M2,M3;
38: PetscScalar M4,m3;
39: } FSubctx;
43: static PetscErrorCode MatFSMult(Mat M ,Vec x,Vec y)
44: {
46: FSubctx *ctx;
47: PetscScalar t;
48:
50: MatShellGetContext(M,&ctx);
51: VecDot(x,ctx->M3,&t);
52: t *= ctx->m3/ctx->M4;
53: MatMult(ctx->M1,x,y);
54: VecAXPY(y,-t,ctx->M2);
55: return(0);
56: }
60: static PetscErrorCode PEPSimpleNRefSetUp(PEP pep,PEPSimpNRefctx **ctx_)
61: {
63: PetscInt i,si,j,n0,m0,nloc,*idx1,*idx2,ne;
64: IS is1,is2;
65: PEPSimpNRefctx *ctx;
66: Vec v;
67: PetscMPIInt rank,size;
70: PetscCalloc1(1,ctx_);
71: ctx = *ctx_;
72: if (pep->npart==1) {
73: pep->refinesubc = NULL;
74: ctx->scatter_id = NULL;
75: ctx->A = pep->A;
76: } else {
77: PetscMalloc2(pep->nmat,&ctx->A,pep->npart,&ctx->scatter_id);
79: /* Duplicate matrices */
80: for (i=0;i<pep->nmat;i++) {
81: MatCreateRedundantMatrix(pep->A[i],0,PetscSubcommChild(pep->refinesubc),MAT_INITIAL_MATRIX,&ctx->A[i]);
82: }
83: MatCreateVecs(ctx->A[0],&ctx->v,NULL);
85: /* Create scatters for sending vectors to each subcommucator */
86: BVGetColumn(pep->V,0,&v);
87: VecGetOwnershipRange(v,&n0,&m0);
88: BVRestoreColumn(pep->V,0,&v);
89: VecGetLocalSize(ctx->v,&nloc);
90: PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
91: VecCreateMPI(PetscObjectComm((PetscObject)pep),nloc,PETSC_DECIDE,&ctx->vg);
92: for (si=0;si<pep->npart;si++) {
93: j = 0;
94: for (i=n0;i<m0;i++) {
95: idx1[j] = i;
96: idx2[j++] = i+pep->n*si;
97: }
98: ISCreateGeneral(PetscObjectComm((PetscObject)pep),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
99: ISCreateGeneral(PetscObjectComm((PetscObject)pep),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
100: BVGetColumn(pep->V,0,&v);
101: VecScatterCreate(v,is1,ctx->vg,is2,&ctx->scatter_id[si]);
102: BVRestoreColumn(pep->V,0,&v);
103: ISDestroy(&is1);
104: ISDestroy(&is2);
105: }
106: PetscFree2(idx1,idx2);
107: }
108: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT){
109: MPI_Comm_rank(PetscObjectComm((PetscObject)ctx->A[0]),&rank);
110: MPI_Comm_size(PetscObjectComm((PetscObject)ctx->A[0]),&size);
111: if (size>1) {
112: if (pep->npart==1) {
113: BVGetColumn(pep->V,0,&v);
114: } else v = ctx->v;
115: VecGetOwnershipRange(v,&n0,&m0);
116: ne = (rank == size-1)?pep->n:0;
117: VecCreateMPI(PetscObjectComm((PetscObject)ctx->A[0]),ne,PETSC_DECIDE,&ctx->nv);
118: PetscMalloc1(m0-n0,&idx1);
119: for (i=n0;i<m0;i++) idx1[i-n0] = i;
120: ISCreateGeneral(PetscObjectComm((PetscObject)pep),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
121: VecScatterCreate(v,is1,ctx->nv,is1,&ctx->nst);
122: if (pep->npart==1) {
123: BVRestoreColumn(pep->V,0,&v);
124: }
125: PetscFree(idx1);
126: ISDestroy(&is1);
127: }
128: }
129: return(0);
130: }
132: /*
133: Gather Eigenpair idx from subcommunicator with color sc
134: */
137: static PetscErrorCode PEPSimpleNRefGatherEigenpair(PEP pep,PEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx,PetscInt *fail)
138: {
139: PetscErrorCode ierr;
140: PetscMPIInt nproc,p;
141: MPI_Comm comm=((PetscObject)pep)->comm;
142: Vec v;
143: const PetscScalar *array;
146: if (pep->npart>1) {
147: MPI_Comm_size(comm,&nproc);
148: p = (nproc/pep->npart)*sc+PetscMin(sc,nproc%pep->npart);
149: /* Communicate convergence successful */
150: MPI_Bcast(fail,1,MPIU_INT,p,comm);
151: if (!(*fail)) {
152: /* Process 0 of subcommunicator sc broadcasts the eigenvalue */
153: MPI_Bcast(&pep->eigr[idx],1,MPIU_SCALAR,p,comm);
154: /* Gather pep->V[idx] from the subcommuniator sc */
155: BVGetColumn(pep->V,idx,&v);
156: if (pep->refinesubc->color==sc) {
157: VecGetArrayRead(ctx->v,&array);
158: VecPlaceArray(ctx->vg,array);
159: }
160: VecScatterBegin(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
161: VecScatterEnd(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
162: if (pep->refinesubc->color==sc) {
163: VecResetArray(ctx->vg);
164: VecRestoreArrayRead(ctx->v,&array);
165: }
166: BVRestoreColumn(pep->V,idx,&v);
167: }
168: } else {
169: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT && !(*fail)) {
170: MPI_Comm_size(comm,&nproc);
171: p = (nproc/pep->npart)*sc+PetscMin(sc,nproc%pep->npart);
172: MPI_Bcast(&pep->eigr[idx],1,MPIU_SCALAR,p,comm);
173: }
174: }
175: return(0);
176: }
180: static PetscErrorCode PEPSimpleNRefScatterEigenvector(PEP pep,PEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx)
181: {
182: PetscErrorCode ierr;
183: Vec v;
184: const PetscScalar *array;
185:
187: if (pep->npart>1) {
188: BVGetColumn(pep->V,idx,&v);
189: if (pep->refinesubc->color==sc) {
190: VecGetArrayRead(ctx->v,&array);
191: VecPlaceArray(ctx->vg,array);
192: }
193: VecScatterBegin(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
194: VecScatterEnd(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
195: if (pep->refinesubc->color==sc) {
196: VecResetArray(ctx->vg);
197: VecRestoreArrayRead(ctx->v,&array);
198: }
199: BVRestoreColumn(pep->V,idx,&v);
200: }
201: return(0);
202: }
206: static PetscErrorCode PEPEvaluateFunctionDerivatives(PEP pep,PetscScalar alpha,PetscScalar *vals)
207: {
208: PetscInt i,nmat=pep->nmat;
209: PetscScalar a0,a1,a2;
210: PetscReal *a=pep->pbc,*b=a+nmat,*g=b+nmat;
213: a0 = 0.0;
214: a1 = 1.0;
215: vals[0] = 0.0;
216: if (nmat>1) vals[1] = 1/a[0];
217: for (i=2;i<nmat;i++) {
218: a2 = ((alpha-b[i-2])*a1-g[i-2]*a0)/a[i-2];
219: vals[i] = (a2+(alpha-b[i-1])*vals[i-1]-g[i-1]*vals[i-2])/a[i-1];
220: a0 = a1; a1 = a2;
221: }
222: return(0);
223: }
227: static PetscErrorCode PEPSimpleNRefSetUpSystem(PEP pep,Mat *A,PEPSimpNRefctx *ctx,PetscInt idx,Mat *Mt,Mat *T,Mat *P,PetscBool ini,Vec t,Vec v)
228: {
229: PetscErrorCode ierr;
230: PetscInt i,nmat=pep->nmat,ml,m0,n0,m1,mg;
231: PetscInt *dnz,*onz,ncols,*cols2=NULL,*nnz;
232: PetscScalar zero=0.0,*coeffs,*coeffs2;
233: PetscMPIInt rank,size;
234: MPI_Comm comm;
235: const PetscInt *cols;
236: const PetscScalar *vals,*array;
237: MatStructure str;
238: FSubctx *fctx;
239: Vec w=ctx->w;
240: Mat M;
243: STGetMatStructure(pep->st,&str);
244: PetscMalloc2(nmat,&coeffs,nmat,&coeffs2);
245: switch (pep->scheme) {
246: case PEP_REFINE_SCHEME_SCHUR:
247: if (ini) {
248: PetscCalloc1(1,&fctx);
249: MatGetSize(A[0],&m0,&n0);
250: MatCreateShell(PetscObjectComm((PetscObject)A[0]),PETSC_DECIDE,PETSC_DECIDE,m0,n0,fctx,T);
251: MatShellSetOperation(*T,MATOP_MULT,(void(*)(void))MatFSMult);
252: } else {
253: MatShellGetContext(*T,&fctx);
254: }
255: M=fctx->M1;
256: break;
257: case PEP_REFINE_SCHEME_MBE:
258: M=*T;
259: break;
260: case PEP_REFINE_SCHEME_EXPLICIT:
261: M=*Mt;
262: break;
263: }
264: if (ini) {
265: MatDuplicate(A[0],MAT_COPY_VALUES,&M);
266: } else {
267: MatCopy(A[0],M,DIFFERENT_NONZERO_PATTERN);
268: }
269: PEPEvaluateBasis(pep,pep->eigr[idx],0,coeffs,NULL);
270: MatScale(M,coeffs[0]);
271: for (i=1;i<nmat;i++) {
272: MatAXPY(M,coeffs[i],A[i],(ini)?str:SUBSET_NONZERO_PATTERN);
273: }
274: PEPEvaluateFunctionDerivatives(pep,pep->eigr[idx],coeffs2);
275: for (i=0;i<nmat && PetscAbsScalar(coeffs2[i])==0.0;i++);
276: MatMult(A[i],v,w);
277: if (coeffs2[i]!=1.0) {
278: VecScale(w,coeffs2[i]);
279: }
280: for (i++;i<nmat;i++) {
281: MatMult(A[i],v,t);
282: VecAXPY(w,coeffs2[i],t);
283: }
284: switch (pep->scheme) {
285: case PEP_REFINE_SCHEME_EXPLICIT:
286: comm = PetscObjectComm((PetscObject)A[0]);
287: MPI_Comm_rank(comm,&rank);
288: MPI_Comm_size(comm,&size);
289: MatGetSize(M,&mg,NULL);
290: MatGetOwnershipRange(M,&m0,&m1);
291: if (ini) {
292: MatCreate(comm,T);
293: MatGetLocalSize(M,&ml,NULL);
294: if (rank==size-1) ml++;
295: MatSetSizes(*T,ml,ml,mg+1,mg+1);
296: MatSetFromOptions(*T);
297: MatSetUp(*T);
298: /* Preallocate M */
299: if (size>1) {
300: MatPreallocateInitialize(comm,ml,ml,dnz,onz);
301: for (i=m0;i<m1;i++) {
302: MatGetRow(M,i,&ncols,&cols,NULL);
303: MatPreallocateSet(i,ncols,cols,dnz,onz);
304: MatPreallocateSet(i,1,&mg,dnz,onz);
305: MatRestoreRow(M,i,&ncols,&cols,NULL);
306: }
307: if (rank==size-1) {
308: PetscCalloc1(mg+1,&cols2);
309: for (i=0;i<mg+1;i++) cols2[i]=i;
310: MatPreallocateSet(m1,mg+1,cols2,dnz,onz);
311: PetscFree(cols2);
312: }
313: MatMPIAIJSetPreallocation(*T,0,dnz,0,onz);
314: MatPreallocateFinalize(dnz,onz);
315: } else {
316: PetscCalloc1(mg+1,&nnz);
317: for (i=0;i<mg;i++) {
318: MatGetRow(M,i,&ncols,NULL,NULL);
319: nnz[i] = ncols+1;
320: MatRestoreRow(M,i,&ncols,NULL,NULL);
321: }
322: nnz[mg] = mg+1;
323: MatSeqAIJSetPreallocation(*T,0,nnz);
324: PetscFree(nnz);
325: }
326: *Mt = M;
327: *P = *T;
328: }
329:
330: /* Set values */
331: VecGetArrayRead(w,&array);
332: for (i=m0;i<m1;i++) {
333: MatGetRow(M,i,&ncols,&cols,&vals);
334: MatSetValues(*T,1,&i,ncols,cols,vals,INSERT_VALUES);
335: MatRestoreRow(M,i,&ncols,&cols,&vals);
336: MatSetValues(*T,1,&i,1,&mg,array+i-m0,INSERT_VALUES);
337: }
338: VecRestoreArrayRead(w,&array);
339: VecConjugate(v);
340: MPI_Comm_size(PetscObjectComm((PetscObject)A[0]),&size);
341: MPI_Comm_rank(PetscObjectComm((PetscObject)A[0]),&rank);
342: if (size>1) {
343: if (rank==size-1) {
344: PetscMalloc1(pep->n,&cols2);
345: for (i=0;i<pep->n;i++) cols2[i]=i;
346: }
347: VecScatterBegin(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD);
348: VecScatterEnd(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD);
349: VecGetArrayRead(ctx->nv,&array);
350: if (rank==size-1) {
351: MatSetValues(*T,1,&mg,pep->n,cols2,array,INSERT_VALUES);
352: MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES);
353: }
354: VecRestoreArrayRead(ctx->nv,&array);
355: } else {
356: PetscMalloc1(m1-m0,&cols2);
357: for (i=0;i<m1-m0;i++) cols2[i]=m0+i;
358: VecGetArrayRead(v,&array);
359: MatSetValues(*T,1,&mg,m1-m0,cols2,array,INSERT_VALUES);
360: MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES);
361: VecRestoreArrayRead(v,&array);
362: }
363: VecConjugate(v);
364: MatAssemblyBegin(*T,MAT_FINAL_ASSEMBLY);
365: MatAssemblyEnd(*T,MAT_FINAL_ASSEMBLY);
366: PetscFree(cols2);
367: break;
368: case PEP_REFINE_SCHEME_SCHUR:
369: fctx->M2 = ctx->w;
370: fctx->M3 = v;
371: fctx->m3 = 0.0;
372: for (i=1;i<nmat-1;i++) fctx->m3 += PetscConj(coeffs[i])*coeffs[i];
373: fctx->M4 = 0.0;
374: for (i=1;i<nmat-1;i++) fctx->M4 += PetscConj(coeffs[i])*coeffs2[i];
375: fctx->M1 = M;
376: if (ini) {
377: MatDuplicate(M,MAT_COPY_VALUES,P);
378: } else {
379: MatCopy(M,*P,SAME_NONZERO_PATTERN);
380: }
381: VecConjugate(v);
382: VecPointwiseMult(t,v,w);
383: VecConjugate(v);
384: VecScale(t,-fctx->m3/fctx->M4);
385: MatDiagonalSet(*P,t,ADD_VALUES);
386: break;
387: case PEP_REFINE_SCHEME_MBE:
388: *T = M;
389: *P = M;
390: break;
391: }
392: PetscFree2(coeffs,coeffs2);
393: return(0);
394: }
398: PetscErrorCode PEPNewtonRefinementSimple(PEP pep,PetscInt *maxits,PetscReal tol,PetscInt k)
399: {
400: PetscErrorCode ierr;
401: PetscInt i,n,its,idx=0,*idx_sc,*its_sc,color,*fail_sc;
402: PetscMPIInt rank,size;
403: Mat Mt=NULL,T=NULL,P=NULL;
404: MPI_Comm comm;
405: Vec r,v,dv,rr=NULL,dvv=NULL,t[2];
406: PetscScalar *array2,deig=0.0,tt[2],ttt;
407: const PetscScalar *array;
408: PetscReal norm,error;
409: PetscBool ini=PETSC_TRUE,sc_pend,solved=PETSC_FALSE;
410: PEPSimpNRefctx *ctx;
411: FSubctx *fctx=NULL;
412: KSPConvergedReason reason;
415: PetscLogEventBegin(PEP_Refine,pep,0,0,0);
416: PEPSimpleNRefSetUp(pep,&ctx);
417: its = (maxits)?*maxits:NREF_MAXIT;
418: if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
419: comm = (pep->npart==1)?PetscObjectComm((PetscObject)pep):PetscSubcommChild(pep->refinesubc);
420: if (pep->npart==1) {
421: BVGetColumn(pep->V,0,&v);
422: } else v = ctx->v;
423: VecDuplicate(v,&ctx->w);
424: VecDuplicate(v,&r);
425: VecDuplicate(v,&dv);
426: VecDuplicate(v,&t[0]);
427: VecDuplicate(v,&t[1]);
428: if (pep->npart==1) { BVRestoreColumn(pep->V,0,&v); }
429: MPI_Comm_size(comm,&size);
430: MPI_Comm_rank(comm,&rank);
431: VecGetLocalSize(r,&n);
432: PetscMalloc3(pep->npart,&idx_sc,pep->npart,&its_sc,pep->npart,&fail_sc);
433: for (i=0;i<pep->npart;i++) fail_sc[i] = 0;
434: for (i=0;i<pep->npart;i++) its_sc[i] = 0;
435: color = (pep->npart==1)?0:pep->refinesubc->color;
436:
437: /* Loop performing iterative refinements */
438: while (!solved) {
439: for (i=0;i<pep->npart;i++) {
440: sc_pend = PETSC_TRUE;
441: if (its_sc[i]==0) {
442: idx_sc[i] = idx++;
443: if (idx_sc[i]>=k) {
444: sc_pend = PETSC_FALSE;
445: } else {
446: PEPSimpleNRefScatterEigenvector(pep,ctx,i,idx_sc[i]);
447: }
448: } else { /* Gather Eigenpair from subcommunicator i */
449: PEPSimpleNRefGatherEigenpair(pep,ctx,i,idx_sc[i],&fail_sc[i]);
450: }
451: while (sc_pend) {
452: if (!fail_sc[i]) {
453: PEPComputeError(pep,idx_sc[i],PEP_ERROR_BACKWARD,&error);
454: }
455: if (error<=tol || its_sc[i]>=its || fail_sc[i]) {
456: idx_sc[i] = idx++;
457: its_sc[i] = 0;
458: if (idx_sc[i]<k) { PEPSimpleNRefScatterEigenvector(pep,ctx,i,idx_sc[i]); }
459: } else {
460: sc_pend = PETSC_FALSE;
461: its_sc[i]++;
462: }
463: if (idx_sc[i]>=k) sc_pend = PETSC_FALSE;
464: }
465: }
466: solved = PETSC_TRUE;
467: for (i=0;i<pep->npart&&solved;i++) solved = PetscNot(idx_sc[i]<k);
468: if (idx_sc[color]<k) {
469: #if !defined(PETSC_USE_COMPLEX)
470: if (pep->eigi[idx_sc[color]]!=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Simple Refinement not implemented in real scalars for complex eigenvalues");
471: #endif
472: if (pep->npart==1) {
473: BVGetColumn(pep->V,idx_sc[color],&v);
474: } else v = ctx->v;
475: PEPSimpleNRefSetUpSystem(pep,ctx->A,ctx,idx_sc[color],&Mt,&T,&P,ini,t[0],v);
476: KSPSetOperators(pep->refineksp,T,P);
477: if (ini) {
478: KSPSetFromOptions(pep->refineksp);
479: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT) {
480: MatCreateVecs(T,&dvv,NULL);
481: VecDuplicate(dvv,&rr);
482: }
483: ini = PETSC_FALSE;
484: }
486: switch (pep->scheme) {
487: case PEP_REFINE_SCHEME_EXPLICIT:
488: MatMult(Mt,v,r);
489: VecGetArrayRead(r,&array);
490: if (rank==size-1) {
491: VecGetArray(rr,&array2);
492: PetscMemcpy(array2,array,n*sizeof(PetscScalar));
493: array2[n] = 0.0;
494: VecRestoreArray(rr,&array2);
495: } else {
496: VecPlaceArray(rr,array);
497: }
498: KSPSolve(pep->refineksp,rr,dvv);
499: KSPGetConvergedReason(pep->refineksp,&reason);
500: if (reason>0) {
501: if (rank != size-1) {
502: VecResetArray(rr);
503: }
504: VecRestoreArrayRead(r,&array);
505: VecGetArrayRead(dvv,&array);
506: VecPlaceArray(dv,array);
507: VecAXPY(v,-1.0,dv);
508: VecNorm(v,NORM_2,&norm);
509: VecScale(v,1.0/norm);
510: VecResetArray(dv);
511: if (rank==size-1) pep->eigr[idx_sc[color]] -= array[n];
512: VecRestoreArrayRead(dvv,&array);
513: } else fail_sc[color] = 1;
514: break;
515: case PEP_REFINE_SCHEME_MBE:
516: MatMult(T,v,r);
517: /* Mixed block elimination */
518: VecConjugate(v);
519: KSPSolveTranspose(pep->refineksp,v,t[0]);
520: KSPGetConvergedReason(pep->refineksp,&reason);
521: if (reason>0) {
522: VecConjugate(t[0]);
523: VecDot(ctx->w,t[0],&tt[0]);
524: KSPSolve(pep->refineksp,ctx->w,t[1]);
525: KSPGetConvergedReason(pep->refineksp,&reason);
526: if (reason>0) {
527: VecDot(t[1],v,&tt[1]);
528: VecDot(r,t[0],&ttt);
529: tt[0] = ttt/tt[0];
530: VecAXPY(r,-tt[0],ctx->w);
531: KSPSolve(pep->refineksp,r,dv);
532: KSPGetConvergedReason(pep->refineksp,&reason);
533: if (reason>0) {
534: VecDot(dv,v,&ttt);
535: tt[1] = ttt/tt[1];
536: VecAXPY(dv,-tt[1],t[1]);
537: deig = tt[0]+tt[1];
538: }
539: }
540: VecConjugate(v);
541: VecAXPY(v,-1.0,dv);
542: VecNorm(v,NORM_2,&norm);
543: VecScale(v,1.0/norm);
544: pep->eigr[idx_sc[color]] -= deig;
545: fail_sc[color] = 0;
546: } else {
547: VecConjugate(v);
548: fail_sc[color] = 1;
549: }
550: break;
551: case PEP_REFINE_SCHEME_SCHUR:
552: MatShellGetContext(T,&fctx);
553: MatMult(fctx->M1,v,r);
554: KSPSolve(pep->refineksp,r,dv);
555: KSPGetConvergedReason(pep->refineksp,&reason);
556: if (reason>0) {
557: VecDot(dv,v,&deig);
558: deig *= -fctx->m3/fctx->M4;
559: VecAXPY(v,-1.0,dv);
560: VecNorm(v,NORM_2,&norm);
561: VecScale(v,1.0/norm);
562: pep->eigr[idx_sc[color]] -= deig;
563: fail_sc[color] = 0;
564: } else fail_sc[color] = 1;
565: break;
566: }
567: if (pep->npart==1) { BVRestoreColumn(pep->V,idx_sc[color],&v); }
568: }
569: }
570: VecDestroy(&t[0]);
571: VecDestroy(&t[1]);
572: VecDestroy(&dv);
573: VecDestroy(&ctx->w);
574: VecDestroy(&r);
575: PetscFree3(idx_sc,its_sc,fail_sc);
576: VecScatterDestroy(&ctx->nst);
577: if (pep->npart>1) {
578: VecDestroy(&ctx->vg);
579: VecDestroy(&ctx->v);
580: for (i=0;i<pep->nmat;i++) {
581: MatDestroy(&ctx->A[i]);
582: }
583: for (i=0;i<pep->npart;i++) {
584: VecScatterDestroy(&ctx->scatter_id[i]);
585: }
586: PetscFree2(ctx->A,ctx->scatter_id);
587: }
588: if (fctx && pep->scheme==PEP_REFINE_SCHEME_SCHUR) {
589: MatDestroy(&P);
590: MatDestroy(&fctx->M1);
591: PetscFree(fctx);
592: }
593: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT) {
594: MatDestroy(&Mt);
595: VecDestroy(&dvv);
596: VecDestroy(&rr);
597: VecDestroy(&ctx->nv);
598: }
599: MatDestroy(&T);
600: PetscFree(ctx);
601: PetscLogEventEnd(PEP_Refine,pep,0,0,0);
602: return(0);
603: }