Actual source code: neprefine.c
slepc-3.18.0 2022-10-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Newton refinement for NEP, simple version
12: */
14: #include <slepc/private/nepimpl.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: FN *fn;
24: } NEPSimpNRefctx;
26: typedef struct {
27: Mat M1;
28: Vec M2,M3;
29: PetscScalar M4,m3;
30: } NEP_REFINE_MATSHELL;
32: static PetscErrorCode MatMult_FS(Mat M ,Vec x,Vec y)
33: {
34: NEP_REFINE_MATSHELL *ctx;
35: PetscScalar t;
37: MatShellGetContext(M,&ctx);
38: VecDot(x,ctx->M3,&t);
39: t *= ctx->m3/ctx->M4;
40: MatMult(ctx->M1,x,y);
41: VecAXPY(y,-t,ctx->M2);
42: return 0;
43: }
45: static PetscErrorCode NEPSimpleNRefSetUp(NEP nep,NEPSimpNRefctx **ctx_)
46: {
47: PetscInt i,si,j,n0,m0,nloc,*idx1,*idx2,ne;
48: IS is1,is2;
49: NEPSimpNRefctx *ctx;
50: Vec v;
51: PetscMPIInt rank,size;
52: MPI_Comm child;
54: PetscCalloc1(1,ctx_);
55: ctx = *ctx_;
56: if (nep->npart==1) {
57: ctx->A = nep->A;
58: ctx->fn = nep->f;
59: nep->refinesubc = NULL;
60: ctx->scatter_id = NULL;
61: } else {
62: PetscSubcommGetChild(nep->refinesubc,&child);
63: PetscMalloc2(nep->nt,&ctx->A,nep->npart,&ctx->scatter_id);
65: /* Duplicate matrices */
66: for (i=0;i<nep->nt;i++) MatCreateRedundantMatrix(nep->A[i],0,child,MAT_INITIAL_MATRIX,&ctx->A[i]);
67: MatCreateVecs(ctx->A[0],&ctx->v,NULL);
69: /* Duplicate FNs */
70: PetscMalloc1(nep->nt,&ctx->fn);
71: for (i=0;i<nep->nt;i++) FNDuplicate(nep->f[i],child,&ctx->fn[i]);
73: /* Create scatters for sending vectors to each subcommucator */
74: BVGetColumn(nep->V,0,&v);
75: VecGetOwnershipRange(v,&n0,&m0);
76: BVRestoreColumn(nep->V,0,&v);
77: VecGetLocalSize(ctx->v,&nloc);
78: PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2);
79: VecCreateMPI(PetscObjectComm((PetscObject)nep),nloc,PETSC_DECIDE,&ctx->vg);
80: for (si=0;si<nep->npart;si++) {
81: j = 0;
82: for (i=n0;i<m0;i++) {
83: idx1[j] = i;
84: idx2[j++] = i+nep->n*si;
85: }
86: ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
87: ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx2,PETSC_COPY_VALUES,&is2);
88: BVGetColumn(nep->V,0,&v);
89: VecScatterCreate(v,is1,ctx->vg,is2,&ctx->scatter_id[si]);
90: BVRestoreColumn(nep->V,0,&v);
91: ISDestroy(&is1);
92: ISDestroy(&is2);
93: }
94: PetscFree2(idx1,idx2);
95: }
96: if (nep->scheme==NEP_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 (nep->npart==1) BVGetColumn(nep->V,0,&v);
101: else v = ctx->v;
102: VecGetOwnershipRange(v,&n0,&m0);
103: ne = (rank == size-1)?nep->n:0;
104: VecCreateMPI(PetscObjectComm((PetscObject)ctx->A[0]),ne,PETSC_DECIDE,&ctx->nv);
105: PetscMalloc1(m0-n0,&idx1);
106: for (i=n0;i<m0;i++) {
107: idx1[i-n0] = i;
108: }
109: ISCreateGeneral(PetscObjectComm((PetscObject)ctx->A[0]),(m0-n0),idx1,PETSC_COPY_VALUES,&is1);
110: VecScatterCreate(v,is1,ctx->nv,is1,&ctx->nst);
111: if (nep->npart==1) BVRestoreColumn(nep->V,0,&v);
112: PetscFree(idx1);
113: ISDestroy(&is1);
114: }
115: } return 0;
116: }
118: /*
119: Gather Eigenpair idx from subcommunicator with color sc
120: */
121: static PetscErrorCode NEPSimpleNRefGatherEigenpair(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx,PetscInt *fail)
122: {
123: PetscMPIInt nproc,p;
124: MPI_Comm comm=((PetscObject)nep)->comm;
125: Vec v;
126: PetscScalar *array;
128: MPI_Comm_size(comm,&nproc);
129: p = (nproc/nep->npart)*(sc+1)+PetscMin(sc+1,nproc%nep->npart)-1;
130: if (nep->npart>1) {
131: /* Communicate convergence successful */
132: MPI_Bcast(fail,1,MPIU_INT,p,comm);
133: if (!(*fail)) {
134: /* Process 0 of subcommunicator sc broadcasts the eigenvalue */
135: MPI_Bcast(&nep->eigr[idx],1,MPIU_SCALAR,p,comm);
136: /* Gather nep->V[idx] from the subcommuniator sc */
137: BVGetColumn(nep->V,idx,&v);
138: if (nep->refinesubc->color==sc) {
139: VecGetArray(ctx->v,&array);
140: VecPlaceArray(ctx->vg,array);
141: }
142: VecScatterBegin(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
143: VecScatterEnd(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE);
144: if (nep->refinesubc->color==sc) {
145: VecResetArray(ctx->vg);
146: VecRestoreArray(ctx->v,&array);
147: }
148: BVRestoreColumn(nep->V,idx,&v);
149: }
150: } else {
151: if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT && !(*fail)) MPI_Bcast(&nep->eigr[idx],1,MPIU_SCALAR,p,comm);
152: }
153: return 0;
154: }
156: static PetscErrorCode NEPSimpleNRefScatterEigenvector(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx)
157: {
158: Vec v;
159: PetscScalar *array;
161: if (nep->npart>1) {
162: BVGetColumn(nep->V,idx,&v);
163: if (nep->refinesubc->color==sc) {
164: VecGetArray(ctx->v,&array);
165: VecPlaceArray(ctx->vg,array);
166: }
167: VecScatterBegin(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
168: VecScatterEnd(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD);
169: if (nep->refinesubc->color==sc) {
170: VecResetArray(ctx->vg);
171: VecRestoreArray(ctx->v,&array);
172: }
173: BVRestoreColumn(nep->V,idx,&v);
174: }
175: return 0;
176: }
178: static PetscErrorCode NEPSimpleNRefSetUpSystem(NEP nep,NEPSimpNRefctx *ctx,Mat *A,PetscInt idx,Mat *Mt,Mat *T,Mat *P,PetscBool ini,Vec t,Vec v)
179: {
180: PetscInt i,st,ml,m0,n0,m1,mg;
181: PetscInt *dnz,*onz,ncols,*cols2=NULL,*nnz,nt=nep->nt;
182: PetscScalar zero=0.0,*coeffs,*coeffs2;
183: PetscMPIInt rank,size;
184: MPI_Comm comm;
185: const PetscInt *cols;
186: const PetscScalar *vals,*array;
187: NEP_REFINE_MATSHELL *fctx;
188: Vec w=ctx->w;
189: Mat M;
191: PetscMalloc2(nt,&coeffs,nt,&coeffs2);
192: switch (nep->scheme) {
193: case NEP_REFINE_SCHEME_SCHUR:
194: if (ini) {
195: PetscCalloc1(1,&fctx);
196: MatGetSize(A[0],&m0,&n0);
197: MatCreateShell(PetscObjectComm((PetscObject)A[0]),PETSC_DECIDE,PETSC_DECIDE,m0,n0,fctx,T);
198: MatShellSetOperation(*T,MATOP_MULT,(void(*)(void))MatMult_FS);
199: } else MatShellGetContext(*T,&fctx);
200: M=fctx->M1;
201: break;
202: case NEP_REFINE_SCHEME_MBE:
203: M=*T;
204: break;
205: case NEP_REFINE_SCHEME_EXPLICIT:
206: M=*Mt;
207: break;
208: }
209: if (ini) MatDuplicate(A[0],MAT_COPY_VALUES,&M);
210: else MatCopy(A[0],M,DIFFERENT_NONZERO_PATTERN);
211: for (i=0;i<nt;i++) FNEvaluateFunction(ctx->fn[i],nep->eigr[idx],coeffs+i);
212: if (coeffs[0]!=1.0) MatScale(M,coeffs[0]);
213: for (i=1;i<nt;i++) MatAXPY(M,coeffs[i],A[i],(ini)?nep->mstr:SUBSET_NONZERO_PATTERN);
214: for (i=0;i<nt;i++) FNEvaluateDerivative(ctx->fn[i],nep->eigr[idx],coeffs2+i);
215: st = 0;
216: for (i=0;i<nt && PetscAbsScalar(coeffs2[i])==0.0;i++) st++;
217: MatMult(A[st],v,w);
218: if (coeffs2[st]!=1.0) VecScale(w,coeffs2[st]);
219: for (i=st+1;i<nt;i++) {
220: MatMult(A[i],v,t);
221: VecAXPY(w,coeffs2[i],t);
222: }
224: switch (nep->scheme) {
225: case NEP_REFINE_SCHEME_EXPLICIT:
226: comm = PetscObjectComm((PetscObject)A[0]);
227: MPI_Comm_rank(comm,&rank);
228: MPI_Comm_size(comm,&size);
229: MatGetSize(M,&mg,NULL);
230: MatGetOwnershipRange(M,&m0,&m1);
231: if (ini) {
232: MatCreate(comm,T);
233: MatGetLocalSize(M,&ml,NULL);
234: if (rank==size-1) ml++;
235: MatSetSizes(*T,ml,ml,mg+1,mg+1);
236: MatSetFromOptions(*T);
237: MatSetUp(*T);
238: /* Preallocate M */
239: if (size>1) {
240: MatPreallocateBegin(comm,ml,ml,dnz,onz);
241: for (i=m0;i<m1;i++) {
242: MatGetRow(M,i,&ncols,&cols,NULL);
243: MatPreallocateSet(i,ncols,cols,dnz,onz);
244: MatPreallocateSet(i,1,&mg,dnz,onz);
245: MatRestoreRow(M,i,&ncols,&cols,NULL);
246: }
247: if (rank==size-1) {
248: PetscCalloc1(mg+1,&cols2);
249: for (i=0;i<mg+1;i++) cols2[i]=i;
250: MatPreallocateSet(m1,mg+1,cols2,dnz,onz);
251: PetscFree(cols2);
252: }
253: MatMPIAIJSetPreallocation(*T,0,dnz,0,onz);
254: MatPreallocateEnd(dnz,onz);
255: } else {
256: PetscCalloc1(mg+1,&nnz);
257: for (i=0;i<mg;i++) {
258: MatGetRow(M,i,&ncols,NULL,NULL);
259: nnz[i] = ncols+1;
260: MatRestoreRow(M,i,&ncols,NULL,NULL);
261: }
262: nnz[mg] = mg+1;
263: MatSeqAIJSetPreallocation(*T,0,nnz);
264: PetscFree(nnz);
265: }
266: *Mt = M;
267: *P = *T;
268: }
270: /* Set values */
271: VecGetArrayRead(w,&array);
272: for (i=m0;i<m1;i++) {
273: MatGetRow(M,i,&ncols,&cols,&vals);
274: MatSetValues(*T,1,&i,ncols,cols,vals,INSERT_VALUES);
275: MatRestoreRow(M,i,&ncols,&cols,&vals);
276: MatSetValues(*T,1,&i,1,&mg,array+i-m0,INSERT_VALUES);
277: }
278: VecRestoreArrayRead(w,&array);
279: VecConjugate(v);
280: MPI_Comm_size(PetscObjectComm((PetscObject)A[0]),&size);
281: MPI_Comm_rank(PetscObjectComm((PetscObject)A[0]),&rank);
282: if (size>1) {
283: if (rank==size-1) {
284: PetscMalloc1(nep->n,&cols2);
285: for (i=0;i<nep->n;i++) cols2[i]=i;
286: }
287: VecScatterBegin(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD);
288: VecScatterEnd(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD);
289: VecGetArrayRead(ctx->nv,&array);
290: if (rank==size-1) {
291: MatSetValues(*T,1,&mg,nep->n,cols2,array,INSERT_VALUES);
292: MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES);
293: }
294: VecRestoreArrayRead(ctx->nv,&array);
295: } else {
296: PetscMalloc1(m1-m0,&cols2);
297: for (i=0;i<m1-m0;i++) cols2[i]=m0+i;
298: VecGetArrayRead(v,&array);
299: MatSetValues(*T,1,&mg,m1-m0,cols2,array,INSERT_VALUES);
300: MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES);
301: VecRestoreArrayRead(v,&array);
302: }
303: VecConjugate(v);
304: MatAssemblyBegin(*T,MAT_FINAL_ASSEMBLY);
305: MatAssemblyEnd(*T,MAT_FINAL_ASSEMBLY);
306: PetscFree(cols2);
307: break;
308: case NEP_REFINE_SCHEME_SCHUR:
309: fctx->M2 = ctx->w;
310: fctx->M3 = v;
311: fctx->m3 = 1.0+PetscConj(nep->eigr[idx])*nep->eigr[idx];
312: fctx->M4 = PetscConj(nep->eigr[idx]);
313: fctx->M1 = M;
314: if (ini) MatDuplicate(M,MAT_COPY_VALUES,P);
315: else MatCopy(M,*P,SAME_NONZERO_PATTERN);
316: if (fctx->M4!=0.0) {
317: VecConjugate(v);
318: VecPointwiseMult(t,v,w);
319: VecConjugate(v);
320: VecScale(t,-fctx->m3/fctx->M4);
321: MatDiagonalSet(*P,t,ADD_VALUES);
322: }
323: break;
324: case NEP_REFINE_SCHEME_MBE:
325: *T = M;
326: *P = M;
327: break;
328: }
329: PetscFree2(coeffs,coeffs2);
330: return 0;
331: }
333: PetscErrorCode NEPNewtonRefinementSimple(NEP nep,PetscInt *maxits,PetscReal tol,PetscInt k)
334: {
335: PetscInt i,n,its,idx=0,*idx_sc,*its_sc,color,*fail_sc;
336: PetscMPIInt rank,size;
337: Mat Mt=NULL,T=NULL,P=NULL;
338: MPI_Comm comm;
339: Vec r,v,dv,rr=NULL,dvv=NULL,t[2];
340: const PetscScalar *array;
341: PetscScalar *array2,deig=0.0,tt[2],ttt;
342: PetscReal norm,error;
343: PetscBool ini=PETSC_TRUE,sc_pend,solved=PETSC_FALSE;
344: NEPSimpNRefctx *ctx;
345: NEP_REFINE_MATSHELL *fctx=NULL;
346: KSPConvergedReason reason;
348: PetscLogEventBegin(NEP_Refine,nep,0,0,0);
349: NEPSimpleNRefSetUp(nep,&ctx);
350: its = (maxits)?*maxits:NREF_MAXIT;
351: if (!nep->refineksp) NEPRefineGetKSP(nep,&nep->refineksp);
352: if (nep->npart==1) BVGetColumn(nep->V,0,&v);
353: else v = ctx->v;
354: VecDuplicate(v,&ctx->w);
355: VecDuplicate(v,&r);
356: VecDuplicate(v,&dv);
357: VecDuplicate(v,&t[0]);
358: VecDuplicate(v,&t[1]);
359: if (nep->npart==1) {
360: BVRestoreColumn(nep->V,0,&v);
361: PetscObjectGetComm((PetscObject)nep,&comm);
362: } else PetscSubcommGetChild(nep->refinesubc,&comm);
363: MPI_Comm_size(comm,&size);
364: MPI_Comm_rank(comm,&rank);
365: VecGetLocalSize(r,&n);
366: PetscMalloc3(nep->npart,&idx_sc,nep->npart,&its_sc,nep->npart,&fail_sc);
367: for (i=0;i<nep->npart;i++) fail_sc[i] = 0;
368: for (i=0;i<nep->npart;i++) its_sc[i] = 0;
369: color = (nep->npart==1)?0:nep->refinesubc->color;
371: /* Loop performing iterative refinements */
372: while (!solved) {
373: for (i=0;i<nep->npart;i++) {
374: sc_pend = PETSC_TRUE;
375: if (its_sc[i]==0) {
376: idx_sc[i] = idx++;
377: if (idx_sc[i]>=k) {
378: sc_pend = PETSC_FALSE;
379: } else NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i]);
380: } else { /* Gather Eigenpair from subcommunicator i */
381: NEPSimpleNRefGatherEigenpair(nep,ctx,i,idx_sc[i],&fail_sc[i]);
382: }
383: while (sc_pend) {
384: if (!fail_sc[i]) NEPComputeError(nep,idx_sc[i],NEP_ERROR_RELATIVE,&error);
385: if (error<=tol || its_sc[i]>=its || fail_sc[i]) {
386: idx_sc[i] = idx++;
387: its_sc[i] = 0;
388: fail_sc[i] = 0;
389: if (idx_sc[i]<k) NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i]);
390: } else {
391: sc_pend = PETSC_FALSE;
392: its_sc[i]++;
393: }
394: if (idx_sc[i]>=k) sc_pend = PETSC_FALSE;
395: }
396: }
397: solved = PETSC_TRUE;
398: for (i=0;i<nep->npart&&solved;i++) solved = PetscNot(idx_sc[i]<k);
399: if (idx_sc[color]<k) {
400: #if !defined(PETSC_USE_COMPLEX)
402: #endif
403: if (nep->npart==1) BVGetColumn(nep->V,idx_sc[color],&v);
404: else v = ctx->v;
405: NEPSimpleNRefSetUpSystem(nep,ctx,ctx->A,idx_sc[color],&Mt,&T,&P,ini,t[0],v);
406: NEP_KSPSetOperators(nep->refineksp,T,P);
407: if (ini) {
408: KSPSetFromOptions(nep->refineksp);
409: if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT) {
410: MatCreateVecs(T,&dvv,NULL);
411: VecDuplicate(dvv,&rr);
412: }
413: ini = PETSC_FALSE;
414: }
415: switch (nep->scheme) {
416: case NEP_REFINE_SCHEME_EXPLICIT:
417: MatMult(Mt,v,r);
418: VecGetArrayRead(r,&array);
419: if (rank==size-1) {
420: VecGetArray(rr,&array2);
421: PetscArraycpy(array2,array,n);
422: array2[n] = 0.0;
423: VecRestoreArray(rr,&array2);
424: } else VecPlaceArray(rr,array);
425: KSPSolve(nep->refineksp,rr,dvv);
426: KSPGetConvergedReason(nep->refineksp,&reason);
427: if (reason>0) {
428: if (rank != size-1) VecResetArray(rr);
429: VecRestoreArrayRead(r,&array);
430: VecGetArrayRead(dvv,&array);
431: VecPlaceArray(dv,array);
432: VecAXPY(v,-1.0,dv);
433: VecNorm(v,NORM_2,&norm);
434: VecScale(v,1.0/norm);
435: VecResetArray(dv);
436: if (rank==size-1) nep->eigr[idx_sc[color]] -= array[n];
437: VecRestoreArrayRead(dvv,&array);
438: } else fail_sc[color] = 1;
439: break;
440: case NEP_REFINE_SCHEME_MBE:
441: MatMult(T,v,r);
442: /* Mixed block elimination */
443: VecConjugate(v);
444: KSPSolveTranspose(nep->refineksp,v,t[0]);
445: KSPGetConvergedReason(nep->refineksp,&reason);
446: if (reason>0) {
447: VecConjugate(t[0]);
448: VecDot(ctx->w,t[0],&tt[0]);
449: KSPSolve(nep->refineksp,ctx->w,t[1]);
450: KSPGetConvergedReason(nep->refineksp,&reason);
451: if (reason>0) {
452: VecDot(t[1],v,&tt[1]);
453: VecDot(r,t[0],&ttt);
454: tt[0] = ttt/tt[0];
455: VecAXPY(r,-tt[0],ctx->w);
456: KSPSolve(nep->refineksp,r,dv);
457: KSPGetConvergedReason(nep->refineksp,&reason);
458: if (reason>0) {
459: VecDot(dv,v,&ttt);
460: tt[1] = ttt/tt[1];
461: VecAXPY(dv,-tt[1],t[1]);
462: deig = tt[0]+tt[1];
463: }
464: }
465: VecConjugate(v);
466: VecAXPY(v,-1.0,dv);
467: VecNorm(v,NORM_2,&norm);
468: VecScale(v,1.0/norm);
469: nep->eigr[idx_sc[color]] -= deig;
470: fail_sc[color] = 0;
471: } else {
472: VecConjugate(v);
473: fail_sc[color] = 1;
474: }
475: break;
476: case NEP_REFINE_SCHEME_SCHUR:
477: fail_sc[color] = 1;
478: MatShellGetContext(T,&fctx);
479: if (fctx->M4!=0.0) {
480: MatMult(fctx->M1,v,r);
481: KSPSolve(nep->refineksp,r,dv);
482: KSPGetConvergedReason(nep->refineksp,&reason);
483: if (reason>0) {
484: VecDot(dv,v,&deig);
485: deig *= -fctx->m3/fctx->M4;
486: VecAXPY(v,-1.0,dv);
487: VecNorm(v,NORM_2,&norm);
488: VecScale(v,1.0/norm);
489: nep->eigr[idx_sc[color]] -= deig;
490: fail_sc[color] = 0;
491: }
492: }
493: break;
494: }
495: if (nep->npart==1) BVRestoreColumn(nep->V,idx_sc[color],&v);
496: }
497: }
498: VecDestroy(&t[0]);
499: VecDestroy(&t[1]);
500: VecDestroy(&dv);
501: VecDestroy(&ctx->w);
502: VecDestroy(&r);
503: PetscFree3(idx_sc,its_sc,fail_sc);
504: VecScatterDestroy(&ctx->nst);
505: if (nep->npart>1) {
506: VecDestroy(&ctx->vg);
507: VecDestroy(&ctx->v);
508: for (i=0;i<nep->nt;i++) MatDestroy(&ctx->A[i]);
509: for (i=0;i<nep->npart;i++) VecScatterDestroy(&ctx->scatter_id[i]);
510: PetscFree2(ctx->A,ctx->scatter_id);
511: }
512: if (fctx && nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
513: MatDestroy(&P);
514: MatDestroy(&fctx->M1);
515: PetscFree(fctx);
516: }
517: if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT) {
518: MatDestroy(&Mt);
519: VecDestroy(&dvv);
520: VecDestroy(&rr);
521: VecDestroy(&ctx->nv);
522: if (nep->npart>1) {
523: for (i=0;i<nep->nt;i++) FNDestroy(&ctx->fn[i]);
524: PetscFree(ctx->fn);
525: }
526: }
527: if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
528: if (nep->npart>1) {
529: for (i=0;i<nep->nt;i++) FNDestroy(&ctx->fn[i]);
530: PetscFree(ctx->fn);
531: }
532: }
533: MatDestroy(&T);
534: PetscFree(ctx);
535: PetscLogEventEnd(NEP_Refine,nep,0,0,0);
536: return 0;
537: }