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