Actual source code: nrefine.c
slepc-3.7.2 2016-07-19
1: /*
2: Newton refinement for polynomial eigenproblems.
4: References:
6: [1] T. Betcke and D. Kressner, "Perturbation, extraction and refinement
7: of invariant pairs for matrix polynomials", Linear Algebra Appl.
8: 435(3):514-536, 2011.
10: [2] C. Campos and J.E. Roman, "Parallel iterative refinement in
11: polynomial eigenvalue problems", submitted, 2015.
13: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14: SLEPc - Scalable Library for Eigenvalue Problem Computations
15: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
17: This file is part of SLEPc.
19: SLEPc is free software: you can redistribute it and/or modify it under the
20: terms of version 3 of the GNU Lesser General Public License as published by
21: the Free Software Foundation.
23: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
24: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
25: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
26: more details.
28: You should have received a copy of the GNU Lesser General Public License
29: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
30: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: */
33: #include <slepc/private/pepimpl.h>
34: #include <slepcblaslapack.h>
36: typedef struct {
37: Mat *A,M1;
38: BV V,M2,M3,W;
39: PetscInt k,nmat;
40: PetscScalar *fih,*work,*M4;
41: PetscBLASInt *pM4;
42: PetscBool compM1;
43: Vec t;
44: } FSubctx;
46: typedef struct {
47: Mat E[2],M1;
48: Vec tN,ttN,t1,vseq;
49: VecScatter scatterctx;
50: PetscBool compM1;
51: PetscInt *map0,*map1,*idxg,*idxp;
52: PetscSubcomm subc;
53: VecScatter scatter_sub;
54: VecScatter *scatter_id,*scatterp_id;
55: Mat *A;
56: BV V,W,M2,M3,Wt;
57: PetscScalar *M4,*w,*wt,*d,*dt;
58: Vec t,tg,Rv,Vi,tp,tpg;
59: PetscInt idx,*cols;
60: } MatExplicitCtx;
64: static PetscErrorCode MatFSMult(Mat M ,Vec x,Vec y)
65: {
66: #if defined(PETSC_MISSING_LAPACK_GETRS)
68: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable");
69: #else
71: FSubctx *ctx;
72: PetscInt k,i;
73: PetscScalar *c;
74: PetscBLASInt k_,one=1,info;
77: MatShellGetContext(M,&ctx);
78: VecCopy(x,ctx->t);
79: k = ctx->k;
80: c = ctx->work;
81: PetscBLASIntCast(k,&k_);
82: MatMult(ctx->M1,x,y);
83: VecConjugate(ctx->t);
84: BVDotVec(ctx->M3,ctx->t,c);
85: for (i=0;i<k;i++) c[i] = PetscConj(c[i]);
86: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&k_,&one,ctx->M4,&k_,ctx->pM4,c,&k_,&info));
87: BVMultVec(ctx->M2,-1.0,1.0,y,c);
88: return(0);
89: #endif
90: }
94: /*
95: Evaluates the first d elements of the polynomial basis
96: on a given matrix H which is considered to be triangular
97: */
98: static PetscErrorCode PEPEvaluateBasisforMatrix(PEP pep,PetscInt nm,PetscInt k,PetscScalar *H,PetscInt ldh,PetscScalar *fH)
99: {
101: PetscInt i,j,ldfh=nm*k,off,nmat=pep->nmat;
102: PetscReal *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat,t;
103: PetscScalar corr=0.0,alpha,beta;
104: PetscBLASInt k_,ldh_,ldfh_;
107: PetscBLASIntCast(ldh,&ldh_);
108: PetscBLASIntCast(k,&k_);
109: PetscBLASIntCast(ldfh,&ldfh_);
110: PetscMemzero(fH,nm*k*k*sizeof(PetscScalar));
111: if (nm>0) for (j=0;j<k;j++) fH[j+j*ldfh] = 1.0;
112: if (nm>1) {
113: t = b[0]/a[0];
114: off = k;
115: for (j=0;j<k;j++) {
116: for (i=0;i<k;i++) fH[off+i+j*ldfh] = H[i+j*ldh]/a[0];
117: fH[j+j*ldfh] -= t;
118: }
119: }
120: for (i=2;i<nm;i++) {
121: off = i*k;
122: if (i==2) {
123: for (j=0;j<k;j++) {
124: fH[off+j+j*ldfh] = 1.0;
125: H[j+j*ldh] -= b[1];
126: }
127: } else {
128: for (j=0;j<k;j++) {
129: PetscMemcpy(fH+off+j*ldfh,fH+(i-2)*k+j*ldfh,k*sizeof(PetscScalar));
130: H[j+j*ldh] += corr-b[i-1];
131: }
132: }
133: corr = b[i-1];
134: beta = -g[i-1]/a[i-1];
135: alpha = 1/a[i-1];
136: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&alpha,H,&ldh_,fH+(i-1)*k,&ldfh_,&beta,fH+off,&ldfh_));
137: }
138: for (j=0;j<k;j++) H[j+j*ldh] += corr;
139: return(0);
140: }
144: static PetscErrorCode NRefSysSetup_shell(PEP pep,PetscInt k,PetscScalar *fH,PetscScalar *S,PetscInt lds,PetscScalar *fh,PetscScalar h,FSubctx *ctx)
145: {
146: #if defined(PETSC_MISSING_LAPACK_GESV)
148: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESV - Lapack routine is unavailable");
149: #else
150: PetscErrorCode ierr;
151: PetscScalar *DHii,*T12,*Tr,*Ts,*array,s,ss,sone=1.0,zero=0.0,*M4=ctx->M4,t,*v,*T;
152: const PetscScalar *m3,*m2;
153: PetscInt i,d,j,nmat=pep->nmat,lda=nmat*k,deg=nmat-1,nloc;
154: PetscReal *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;
155: PetscBLASInt k_,lda_,lds_,nloc_,one=1,info;
156: Mat *A=ctx->A,Mk,M1=ctx->M1,P;
157: BV V=ctx->V,M2=ctx->M2,M3=ctx->M3,W=ctx->W;
158: MatStructure str;
159: Vec vc;
162: STGetMatStructure(pep->st,&str);
163: PetscMalloc3(nmat*k*k,&T12,k*k,&Tr,PetscMax(k*k,nmat),&Ts);
164: DHii = T12;
165: PetscMemzero(DHii,k*k*nmat*sizeof(PetscScalar));
166: for (i=0;i<k;i++) DHii[k+i+i*lda] = 1.0/a[0];
167: for (d=2;d<nmat;d++) {
168: for (j=0;j<k;j++) {
169: for (i=0;i<k;i++) {
170: DHii[d*k+i+j*lda] = ((h-b[d-1])*DHii[(d-1)*k+i+j*lda]+fH[(d-1)*k+i+j*lda]-g[d-1]*DHii[(d-2)*k+i+j*lda])/(a[d-1]);
171: }
172: }
173: }
174: /* T11 */
175: if (!ctx->compM1) {
176: MatCopy(A[0],M1,DIFFERENT_NONZERO_PATTERN);
177: PEPEvaluateBasis(pep,h,0,Ts,NULL);
178: for (j=1;j<nmat;j++) {
179: MatAXPY(M1,Ts[j],A[j],str);
180: }
181: }
183: /* T22 */
184: PetscBLASIntCast(lds,&lds_);
185: PetscBLASIntCast(k,&k_);
186: PetscBLASIntCast(lda,&lda_);
187: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,S,&lds_,S,&lds_,&zero,Tr,&k_));
188: for (i=1;i<deg;i++) {
189: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,Tr,&k_,DHii+i*k,&lda_,&zero,Ts,&k_));
190: s = (i==1)?0.0:1.0;
191: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,fH+i*k,&lda_,Ts,&k_,&s,M4,&k_));
192: }
193: for (i=0;i<k;i++) for (j=0;j<i;j++) { t=M4[i+j*k];M4[i+j*k]=M4[j+i*k];M4[j+i*k]=t; }
195: /* T12 */
196: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&Mk);
197: for (i=1;i<nmat;i++) {
198: MatDenseGetArray(Mk,&array);
199: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,DHii+i*k,&lda_,&zero,array,&k_));
200: MatDenseRestoreArray(Mk,&array);
201: BVSetActiveColumns(W,0,k);
202: BVMult(W,1.0,0.0,V,Mk);
203: if (i==1) {
204: BVMatMult(W,A[i],M2);
205: } else {
206: BVMatMult(W,A[i],M3); /* using M3 as work space */
207: BVMult(M2,1.0,1.0,M3,NULL);
208: }
209: }
211: /* T21 */
212: MatDenseGetArray(Mk,&array);
213: for (i=1;i<deg;i++) {
214: s = (i==1)?0.0:1.0;
215: ss = PetscConj(fh[i]);
216: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&ss,S,&lds_,fH+i*k,&lda_,&s,array,&k_));
217: }
218: MatDenseRestoreArray(Mk,&array);
219: BVSetActiveColumns(M3,0,k);
220: BVMult(M3,1.0,0.0,V,Mk);
221: for (i=0;i<k;i++) {
222: BVGetColumn(M3,i,&vc);
223: VecConjugate(vc);
224: BVRestoreColumn(M3,i,&vc);
225: }
226: MatDestroy(&Mk);
227: PetscFree3(T12,Tr,Ts);
229: VecGetLocalSize(ctx->t,&nloc);
230: PetscBLASIntCast(nloc,&nloc_);
231: PetscMalloc1(nloc*k,&T);
232: KSPGetOperators(pep->refineksp,NULL,&P);
233: if (!ctx->compM1) { MatCopy(ctx->M1,P,SAME_NONZERO_PATTERN); }
234: BVGetArrayRead(ctx->M2,&m2);
235: BVGetArrayRead(ctx->M3,&m3);
236: VecGetArray(ctx->t,&v);
237: for (i=0;i<nloc;i++) for (j=0;j<k;j++) T[j+i*k] = m3[i+j*nloc];
238: PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&nloc_,ctx->M4,&k_,ctx->pM4,T,&k_,&info));
239: for (i=0;i<nloc;i++) v[i] = BLASdot_(&k_,m2+i,&nloc_,T+i*k,&one);
240: VecRestoreArray(ctx->t,&v);
241: BVRestoreArrayRead(ctx->M2,&m2);
242: BVRestoreArrayRead(ctx->M3,&m3);
243: MatDiagonalSet(P,ctx->t,ADD_VALUES);
244: PetscFree(T);
245: KSPSetUp(pep->refineksp);
246: return(0);
247: #endif
248: }
252: static PetscErrorCode NRefSysSolve_shell(KSP ksp,PetscInt nmat,Vec Rv,PetscScalar *Rh,PetscInt k,Vec dVi,PetscScalar *dHi)
253: {
254: #if defined(PETSC_MISSING_LAPACK_GETRS)
256: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable");
257: #else
259: PetscScalar *t0;
260: PetscBLASInt k_,one=1,info,lda_;
261: PetscInt i,lda=nmat*k;
262: Mat M;
263: FSubctx *ctx;
266: KSPGetOperators(ksp,&M,NULL);
267: MatShellGetContext(M,&ctx);
268: PetscCalloc1(k,&t0);
269: PetscBLASIntCast(lda,&lda_);
270: PetscBLASIntCast(k,&k_);
271: for (i=0;i<k;i++) t0[i] = Rh[i];
272: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&k_,&one,ctx->M4,&k_,ctx->pM4,t0,&k_,&info));
273: BVMultVec(ctx->M2,-1.0,1.0,Rv,t0);
274: KSPSolve(ksp,Rv,dVi);
275: VecConjugate(dVi);
276: BVDotVec(ctx->M3,dVi,dHi);
277: VecConjugate(dVi);
278: for (i=0;i<k;i++) dHi[i] = Rh[i]-PetscConj(dHi[i]);
279: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&k_,&one,ctx->M4,&k_,ctx->pM4,dHi,&k_,&info));
280: PetscFree(t0);
281: return(0);
282: #endif
283: }
287: /*
288: Computes the residual P(H,V*S)*e_j for the polynomial
289: */
290: static PetscErrorCode NRefRightSide(PetscInt nmat,PetscReal *pcf,Mat *A,PetscInt k,BV V,PetscScalar *S,PetscInt lds,PetscInt j,PetscScalar *H,PetscInt ldh,PetscScalar *fH,PetscScalar *DfH,PetscScalar *dH,BV dV,PetscScalar *dVS,PetscInt rds,Vec Rv,PetscScalar *Rh,BV W,Vec t)
291: {
293: PetscScalar *DS0,*DS1,*F,beta=0.0,sone=1.0,none=-1.0,tt=0.0,*h,zero=0.0,*Z,*c0;
294: PetscReal *a=pcf,*b=pcf+nmat,*g=b+nmat;
295: PetscInt i,ii,jj,lda;
296: PetscBLASInt lda_,k_,ldh_,lds_,nmat_,k2_,krds_,j_,one=1;
297: Mat M0;
298: Vec w;
301: PetscMalloc4(k*nmat,&h,k*k,&DS0,k*k,&DS1,k*k,&Z);
302: lda = k*nmat;
303: PetscBLASIntCast(k,&k_);
304: PetscBLASIntCast(lds,&lds_);
305: PetscBLASIntCast(lda,&lda_);
306: PetscBLASIntCast(nmat,&nmat_);
307: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&nmat_,&k_,&sone,S,&lds_,fH+j*lda,&k_,&zero,h,&k_));
308: MatCreateSeqDense(PETSC_COMM_SELF,k,nmat,h,&M0);
309: BVSetActiveColumns(W,0,nmat);
310: BVMult(W,1.0,0.0,V,M0);
311: MatDestroy(&M0);
313: BVGetColumn(W,0,&w);
314: MatMult(A[0],w,Rv);
315: BVRestoreColumn(W,0,&w);
316: for (i=1;i<nmat;i++) {
317: BVGetColumn(W,i,&w);
318: MatMult(A[i],w,t);
319: BVRestoreColumn(W,i,&w);
320: VecAXPY(Rv,1.0,t);
321: }
322: /* Update right-hand side */
323: if (j) {
324: PetscBLASIntCast(ldh,&ldh_);
325: PetscMemzero(Z,k*k*sizeof(PetscScalar));
326: PetscMemzero(DS0,k*k*sizeof(PetscScalar));
327: PetscMemcpy(Z+(j-1)*k,dH+(j-1)*k,k*sizeof(PetscScalar));
328: /* Update DfH */
329: for (i=1;i<nmat;i++) {
330: if (i>1) {
331: beta = -g[i-1];
332: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,fH+(i-1)*k,&lda_,Z,&k_,&beta,DS0,&k_));
333: tt += -b[i-1];
334: for (ii=0;ii<k;ii++) H[ii+ii*ldh] += tt;
335: tt = b[i-1];
336: beta = 1.0/a[i-1];
337: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&beta,DS1,&k_,H,&ldh_,&beta,DS0,&k_));
338: F = DS0; DS0 = DS1; DS1 = F;
339: } else {
340: PetscMemzero(DS1,k*k*sizeof(PetscScalar));
341: for (ii=0;ii<k;ii++) DS1[ii+(j-1)*k] = Z[ii+(j-1)*k]/a[0];
342: }
343: for (jj=j;jj<k;jj++) {
344: for (ii=0;ii<k;ii++) DfH[k*i+ii+jj*lda] += DS1[ii+jj*k];
345: }
346: }
347: for (ii=0;ii<k;ii++) H[ii+ii*ldh] += tt;
348: /* Update right-hand side */
349: PetscBLASIntCast(2*k,&k2_);
350: PetscBLASIntCast(j,&j_);
351: PetscBLASIntCast(k+rds,&krds_);
352: c0 = DS0;
353: PetscMemzero(Rh,k*sizeof(PetscScalar));
354: for (i=0;i<nmat;i++) {
355: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&krds_,&j_,&sone,dVS,&k2_,fH+j*lda+i*k,&one,&zero,h,&one));
356: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&k_,&k_,&sone,S,&lds_,DfH+i*k+j*lda,&one,&sone,h,&one));
357: BVMultVec(V,1.0,0.0,t,h);
358: BVSetActiveColumns(dV,0,rds);
359: BVMultVec(dV,1.0,1.0,t,h+k);
360: BVGetColumn(W,i,&w);
361: MatMult(A[i],t,w);
362: BVRestoreColumn(W,i,&w);
363: if (i>0 && i<nmat-1) {
364: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&k_,&k_,&sone,S,&lds_,h,&one,&zero,c0,&one));
365: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&k_,&k_,&none,fH+i*k,&lda_,c0,&one,&sone,Rh,&one));
366: }
367: }
369: for (i=0;i<nmat;i++) h[i] = -1.0;
370: BVMultVec(W,1.0,1.0,Rv,h);
371: }
372: PetscFree4(h,DS0,DS1,Z);
373: return(0);
374: }
378: static PetscErrorCode NRefSysSolve_mbe(PetscInt k,PetscInt sz,BV W,PetscScalar *w,BV Wt,PetscScalar *wt,PetscScalar *d,PetscScalar *dt,KSP ksp,BV T2,BV T3 ,PetscScalar *T4,PetscBool trans,Vec x1,PetscScalar *x2,Vec sol1,PetscScalar *sol2,Vec vw)
379: {
381: PetscInt i,j,incf,incc;
382: PetscScalar *y,*g,*xx2,*ww,y2,*dd;
383: Vec v,t,xx1;
384: BV WW,T;
387: PetscMalloc3(sz,&y,sz,&g,k,&xx2);
388: if (trans) {
389: WW = W; ww = w; dd = d; T = T3; incf = 0; incc = 1;
390: } else {
391: WW = Wt; ww = wt; dd = dt; T = T2; incf = 1; incc = 0;
392: }
393: xx1 = vw;
394: VecCopy(x1,xx1);
395: PetscMemcpy(xx2,x2,sz*sizeof(PetscScalar));
396: PetscMemzero(sol2,k*sizeof(PetscScalar));
397: for (i=sz-1;i>=0;i--) {
398: BVGetColumn(WW,i,&v);
399: VecConjugate(v);
400: VecDot(xx1,v,y+i);
401: VecConjugate(v);
402: BVRestoreColumn(WW,i,&v);
403: for (j=0;j<i;j++) y[i] += ww[j+i*k]*xx2[j];
404: y[i] = -(y[i]-xx2[i])/dd[i];
405: BVGetColumn(T,i,&t);
406: VecAXPY(xx1,-y[i],t);
407: BVRestoreColumn(T,i,&t);
408: for(j=0;j<=i;j++) xx2[j] -= y[i]*T4[j*incf+incc*i+(i*incf+incc*j)*k];
409: g[i] = xx2[i];
410: }
411: if (trans) {
412: KSPSolveTranspose(ksp,xx1,sol1);
413: } else {
414: KSPSolve(ksp,xx1,sol1);
415: }
416: if (trans) {
417: WW = Wt; ww = wt; dd = dt; T = T2; incf = 1; incc = 0;
418: } else {
419: WW = W; ww = w; dd = d; T = T3; incf = 0; incc = 1;
420: }
421: for (i=0;i<sz;i++) {
422: BVGetColumn(T,i,&t);
423: VecConjugate(t);
424: VecDot(sol1,t,&y2);
425: VecConjugate(t);
426: BVRestoreColumn(T,i,&t);
427: for (j=0;j<i;j++) y2 += sol2[j]*T4[j*incf+incc*i+(i*incf+incc*j)*k];
428: y2 = (g[i]-y2)/dd[i];
429: BVGetColumn(WW,i,&v);
430: VecAXPY(sol1,-y2,v);
431: for (j=0;j<i;j++) sol2[j] -= ww[j+i*k]*y2;
432: sol2[i] = y[i]+y2;
433: BVRestoreColumn(WW,i,&v);
434: }
435: PetscFree3(y,g,xx2);
436: return(0);
437: }
441: static PetscErrorCode NRefSysSetup_mbe(PEP pep,PetscInt k,KSP ksp,PetscScalar *fH,PetscScalar *S,PetscInt lds,PetscScalar *fh,PetscScalar h,BV V,MatExplicitCtx *matctx)
442: {
444: PetscInt i,j,l,nmat=pep->nmat,lda=nmat*k,deg=nmat-1;
445: Mat M1=matctx->M1,*A,*At,Mk;
446: PetscReal *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;
447: PetscScalar s,ss,*DHii,*T12,*array,*Ts,*Tr,*M4=matctx->M4,sone=1.0,zero=0.0;
448: PetscScalar *w=matctx->w,*wt=matctx->wt,*d=matctx->d,*dt=matctx->dt;
449: PetscBLASInt lds_,lda_,k_;
450: MatStructure str;
451: PetscBool flg;
452: BV M2=matctx->M2,M3=matctx->M3,W=matctx->W,Wt=matctx->Wt;
453: Vec vc,vc2;
456: PetscMalloc3(nmat*k*k,&T12,k*k,&Tr,PetscMax(k*k,nmat),&Ts);
457: STGetMatStructure(pep->st,&str);
458: STGetTransform(pep->st,&flg);
459: if (flg) {
460: PetscMalloc1(pep->nmat,&At);
461: for (i=0;i<pep->nmat;i++) {
462: STGetTOperators(pep->st,i,&At[i]);
463: }
464: } else At = pep->A;
465: if (matctx->subc) A = matctx->A;
466: else A = At;
467: /* Form the explicit system matrix */
468: DHii = T12;
469: PetscMemzero(DHii,k*k*nmat*sizeof(PetscScalar));
470: for (i=0;i<k;i++) DHii[k+i+i*lda] = 1.0/a[0];
471: for (l=2;l<nmat;l++) {
472: for (j=0;j<k;j++) {
473: for (i=0;i<k;i++) {
474: DHii[l*k+i+j*lda] = ((h-b[l-1])*DHii[(l-1)*k+i+j*lda]+fH[(l-1)*k+i+j*lda]-g[l-1]*DHii[(l-2)*k+i+j*lda])/a[l-1];
475: }
476: }
477: }
479: /* T11 */
480: if (!matctx->compM1) {
481: MatCopy(A[0],M1,DIFFERENT_NONZERO_PATTERN);
482: PEPEvaluateBasis(pep,h,0,Ts,NULL);
483: for (j=1;j<nmat;j++) {
484: MatAXPY(M1,Ts[j],A[j],str);
485: }
486: }
488: /* T22 */
489: PetscBLASIntCast(lds,&lds_);
490: PetscBLASIntCast(k,&k_);
491: PetscBLASIntCast(lda,&lda_);
492: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,S,&lds_,S,&lds_,&zero,Tr,&k_));
493: for (i=1;i<deg;i++) {
494: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,Tr,&k_,DHii+i*k,&lda_,&zero,Ts,&k_));
495: s = (i==1)?0.0:1.0;
496: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,fH+i*k,&lda_,Ts,&k_,&s,M4,&k_));
497: }
499: /* T12 */
500: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&Mk);
501: for (i=1;i<nmat;i++) {
502: MatDenseGetArray(Mk,&array);
503: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,DHii+i*k,&lda_,&zero,array,&k_));
504: MatDenseRestoreArray(Mk,&array);
505: BVSetActiveColumns(W,0,k);
506: BVMult(W,1.0,0.0,V,Mk);
507: if (i==1) {
508: BVMatMult(W,A[i],M2);
509: } else {
510: BVMatMult(W,A[i],M3); /* using M3 as work space */
511: BVMult(M2,1.0,1.0,M3,NULL);
512: }
513: }
515: /* T21 */
516: MatDenseGetArray(Mk,&array);
517: for (i=1;i<deg;i++) {
518: s = (i==1)?0.0:1.0;
519: ss = PetscConj(fh[i]);
520: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&ss,S,&lds_,fH+i*k,&lda_,&s,array,&k_));
521: }
522: MatDenseRestoreArray(Mk,&array);
523: BVSetActiveColumns(M3,0,k);
524: BVMult(M3,1.0,0.0,V,Mk);
525: for (i=0;i<k;i++) {
526: BVGetColumn(M3,i,&vc);
527: VecConjugate(vc);
528: BVRestoreColumn(M3,i,&vc);
529: }
531: KSPSetOperators(ksp,M1,M1);
532: KSPSetUp(ksp);
533: MatDestroy(&Mk);
535: /* Set up for BEMW */
536: for (i=0;i<k;i++) {
537: BVGetColumn(M2,i,&vc);
538: BVGetColumn(W,i,&vc2);
539: NRefSysSolve_mbe(k,i,W,w,Wt,wt,d,dt,ksp,M2,M3,M4,PETSC_FALSE,vc,M4+i*k,vc2,w+i*k,matctx->t);
540: BVRestoreColumn(M2,i,&vc);
541: BVGetColumn(M3,i,&vc);
542: VecConjugate(vc);
543: VecDot(vc2,vc,&d[i]);
544: VecConjugate(vc);
545: BVRestoreColumn(M3,i,&vc);
546: for (j=0;j<i;j++) d[i] += M4[i+j*k]*w[j+i*k];
547: d[i] = M4[i+i*k]-d[i];
548: BVRestoreColumn(W,i,&vc2);
550: BVGetColumn(M3,i,&vc);
551: BVGetColumn(Wt,i,&vc2);
552: for (j=0;j<=i;j++) Ts[j] = M4[i+j*k];
553: NRefSysSolve_mbe(k,i,W,w,Wt,wt,d,dt,ksp,M2,M3,M4,PETSC_TRUE,vc,Ts,vc2,wt+i*k,matctx->t);
554: BVRestoreColumn(M3,i,&vc);
555: BVGetColumn(M2,i,&vc);
556: VecConjugate(vc2);
557: VecDot(vc,vc2,&dt[i]);
558: VecConjugate(vc2);
559: BVRestoreColumn(M2,i,&vc);
560: for (j=0;j<i;j++) dt[i] += M4[j+i*k]*wt[j+i*k];
561: dt[i] = M4[i+i*k]-dt[i];
562: BVRestoreColumn(Wt,i,&vc2);
563: }
565: if (flg) {
566: PetscFree(At);
567: }
568: PetscFree3(T12,Tr,Ts);
569: return(0);
570: }
574: static PetscErrorCode NRefSysSetup_explicit(PEP pep,PetscInt k,KSP ksp,PetscScalar *fH,PetscScalar *S,PetscInt lds,PetscScalar *fh,PetscScalar h,BV V,MatExplicitCtx *matctx,BV W)
575: {
576: PetscErrorCode ierr;
577: PetscInt i,j,d,n,n0,m0,n1,m1,nmat=pep->nmat,lda=nmat*k,deg=nmat-1;
578: PetscInt *idxg=matctx->idxg,*idxp=matctx->idxp,idx,ncols;
579: Mat M,*E=matctx->E,*A,*At,Mk,Md;
580: PetscReal *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;
581: PetscScalar s,ss,*DHii,*T22,*T21,*T12,*Ts,*Tr,*array,*ts,sone=1.0,zero=0.0;
582: PetscBLASInt lds_,lda_,k_;
583: const PetscInt *idxmc;
584: const PetscScalar *valsc,*carray;
585: MatStructure str;
586: Vec vc,vc0;
587: PetscBool flg;
590: PetscMalloc5(k*k,&T22,k*k,&T21,nmat*k*k,&T12,k*k,&Tr,k*k,&Ts);
591: STGetMatStructure(pep->st,&str);
592: KSPGetOperators(ksp,&M,NULL);
593: MatGetOwnershipRange(E[1],&n1,&m1);
594: MatGetOwnershipRange(E[0],&n0,&m0);
595: MatGetOwnershipRange(M,&n,NULL);
596: PetscMalloc1(nmat,&ts);
597: STGetTransform(pep->st,&flg);
598: if (flg) {
599: PetscMalloc1(pep->nmat,&At);
600: for (i=0;i<pep->nmat;i++) {
601: STGetTOperators(pep->st,i,&At[i]);
602: }
603: } else At = pep->A;
604: if (matctx->subc) A = matctx->A;
605: else A = At;
606: /* Form the explicit system matrix */
607: DHii = T12;
608: PetscMemzero(DHii,k*k*nmat*sizeof(PetscScalar));
609: for (i=0;i<k;i++) DHii[k+i+i*lda] = 1.0/a[0];
610: for (d=2;d<nmat;d++) {
611: for (j=0;j<k;j++) {
612: for (i=0;i<k;i++) {
613: DHii[d*k+i+j*lda] = ((h-b[d-1])*DHii[(d-1)*k+i+j*lda]+fH[(d-1)*k+i+j*lda]-g[d-1]*DHii[(d-2)*k+i+j*lda])/a[d-1];
614: }
615: }
616: }
618: /* T11 */
619: if (!matctx->compM1) {
620: MatCopy(A[0],E[0],DIFFERENT_NONZERO_PATTERN);
621: PEPEvaluateBasis(pep,h,0,Ts,NULL);
622: for (j=1;j<nmat;j++) {
623: MatAXPY(E[0],Ts[j],A[j],str);
624: }
625: }
626: for (i=n0;i<m0;i++) {
627: MatGetRow(E[0],i,&ncols,&idxmc,&valsc);
628: idx = n+i-n0;
629: for (j=0;j<ncols;j++) {
630: idxg[j] = matctx->map0[idxmc[j]];
631: }
632: MatSetValues(M,1,&idx,ncols,idxg,valsc,INSERT_VALUES);
633: MatRestoreRow(E[0],i,&ncols,&idxmc,&valsc);
634: }
636: /* T22 */
637: PetscBLASIntCast(lds,&lds_);
638: PetscBLASIntCast(k,&k_);
639: PetscBLASIntCast(lda,&lda_);
640: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,S,&lds_,S,&lds_,&zero,Tr,&k_));
641: for (i=1;i<deg;i++) {
642: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,Tr,&k_,DHii+i*k,&lda_,&zero,Ts,&k_));
643: s = (i==1)?0.0:1.0;
644: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&k_,&sone,fH+i*k,&lda_,Ts,&k_,&s,T22,&k_));
645: }
646: for (j=0;j<k;j++) idxp[j] = matctx->map1[j];
647: for (i=0;i<m1-n1;i++) {
648: idx = n+m0-n0+i;
649: for (j=0;j<k;j++) {
650: Tr[j] = T22[n1+i+j*k];
651: }
652: MatSetValues(M,1,&idx,k,idxp,Tr,INSERT_VALUES);
653: }
655: /* T21 */
656: for (i=1;i<deg;i++) {
657: s = (i==1)?0.0:1.0;
658: ss = PetscConj(fh[i]);
659: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&ss,S,&lds_,fH+i*k,&lda_,&s,T21,&k_));
660: }
661: BVSetActiveColumns(W,0,k);
662: MatCreateSeqDense(PETSC_COMM_SELF,k,k,T21,&Mk);
663: BVMult(W,1.0,0.0,V,Mk);
664: for (i=0;i<k;i++) {
665: BVGetColumn(W,i,&vc);
666: VecConjugate(vc);
667: VecGetArrayRead(vc,&carray);
668: idx = matctx->map1[i];
669: MatSetValues(M,1,&idx,m0-n0,matctx->map0+n0,carray,INSERT_VALUES);
670: VecRestoreArrayRead(vc,&carray);
671: BVRestoreColumn(W,i,&vc);
672: }
674: /* T12 */
675: for (i=1;i<nmat;i++) {
676: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,DHii+i*k,&lda_,&zero,Ts,&k_));
677: for (j=0;j<k;j++) {
678: PetscMemcpy(T12+i*k+j*lda,Ts+j*k,k*sizeof(PetscScalar));
679: }
680: }
681: MatCreateSeqDense(PETSC_COMM_SELF,k,nmat-1,NULL,&Md);
682: for (i=0;i<nmat;i++) ts[i] = 1.0;
683: for (j=0;j<k;j++) {
684: MatDenseGetArray(Md,&array);
685: PetscMemcpy(array,T12+k+j*lda,(nmat-1)*k*sizeof(PetscScalar));
686: MatDenseRestoreArray(Md,&array);
687: BVSetActiveColumns(W,0,nmat-1);
688: BVMult(W,1.0,0.0,V,Md);
689: for (i=nmat-1;i>0;i--) {
690: BVGetColumn(W,i-1,&vc0);
691: BVGetColumn(W,i,&vc);
692: MatMult(A[i],vc0,vc);
693: BVRestoreColumn(W,i-1,&vc0);
694: BVRestoreColumn(W,i,&vc);
695: }
696: BVSetActiveColumns(W,1,nmat);
697: BVGetColumn(W,0,&vc0);
698: BVMultVec(W,1.0,0.0,vc0,ts);
699: VecGetArrayRead(vc0,&carray);
700: idx = matctx->map1[j];
701: MatSetValues(M,m0-n0,matctx->map0+n0,1,&idx,carray,INSERT_VALUES);
702: VecRestoreArrayRead(vc0,&carray);
703: BVRestoreColumn(W,0,&vc0);
704: }
705: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
706: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
707: KSPSetOperators(ksp,M,M);
708: KSPSetUp(ksp);
709: PetscFree(ts);
710: MatDestroy(&Mk);
711: MatDestroy(&Md);
712: if (flg) {
713: PetscFree(At);
714: }
715: PetscFree5(T22,T21,T12,Tr,Ts);
716: return(0);
717: }
721: static PetscErrorCode NRefSysSolve_explicit(PetscInt k,KSP ksp,Vec Rv,PetscScalar *Rh,Vec dVi,PetscScalar *dHi,MatExplicitCtx *matctx)
722: {
723: PetscErrorCode ierr;
724: PetscInt n0,m0,n1,m1,i;
725: PetscScalar *arrayV;
726: const PetscScalar *array;
729: MatGetOwnershipRange(matctx->E[1],&n1,&m1);
730: MatGetOwnershipRange(matctx->E[0],&n0,&m0);
732: /* Right side */
733: VecGetArrayRead(Rv,&array);
734: VecSetValues(matctx->tN,m0-n0,matctx->map0+n0,array,INSERT_VALUES);
735: VecRestoreArrayRead(Rv,&array);
736: VecSetValues(matctx->tN,m1-n1,matctx->map1+n1,Rh+n1,INSERT_VALUES);
737: VecAssemblyBegin(matctx->tN);
738: VecAssemblyEnd(matctx->tN);
740: /* Solve */
741: KSPSolve(ksp,matctx->tN,matctx->ttN);
743: /* Retrieve solution */
744: VecGetArray(dVi,&arrayV);
745: VecGetArrayRead(matctx->ttN,&array);
746: PetscMemcpy(arrayV,array,(m0-n0)*sizeof(PetscScalar));
747: VecRestoreArray(dVi,&arrayV);
748: if (!matctx->subc) {
749: VecGetArray(matctx->t1,&arrayV);
750: for (i=0;i<m1-n1;i++) arrayV[i] = array[m0-n0+i];
751: VecRestoreArray(matctx->t1,&arrayV);
752: VecRestoreArrayRead(matctx->ttN,&array);
753: VecScatterBegin(matctx->scatterctx,matctx->t1,matctx->vseq,INSERT_VALUES,SCATTER_FORWARD);
754: VecScatterEnd(matctx->scatterctx,matctx->t1,matctx->vseq,INSERT_VALUES,SCATTER_FORWARD);
755: VecGetArrayRead(matctx->vseq,&array);
756: for (i=0;i<k;i++) dHi[i] = array[i];
757: VecRestoreArrayRead(matctx->vseq,&array);
758: }
759: return(0);
760: }
764: static PetscErrorCode NRefSysIter(PetscInt i,PEP pep,PetscInt k,KSP ksp,PetscScalar *fH,PetscScalar *S,PetscInt lds,PetscScalar *fh,PetscScalar *H,PetscInt ldh,Vec Rv,PetscScalar *Rh,BV V,Vec dVi,PetscScalar *dHi,MatExplicitCtx *matctx,BV W)
765: {
766: PetscErrorCode ierr;
767: PetscInt j,m,lda=pep->nmat*k,n0,m0,idx;
768: PetscMPIInt root,len;
769: PetscScalar *array2,h;
770: const PetscScalar *array;
771: Vec R,Vi;
772: FSubctx *ctx;
773: Mat M;
776: if (!matctx || !matctx->subc) {
777: for (j=0;j<pep->nmat;j++) fh[j] = fH[j*k+i+i*lda];
778: h = H[i+i*ldh];
779: idx = i;
780: R = Rv;
781: Vi = dVi;
782: switch (pep->scheme) {
783: case PEP_REFINE_SCHEME_EXPLICIT:
784: NRefSysSetup_explicit(pep,k,ksp,fH,S,lds,fh,h,V,matctx,W);
785: matctx->compM1 = PETSC_FALSE;
786: break;
787: case PEP_REFINE_SCHEME_MBE:
788: NRefSysSetup_mbe(pep,k,ksp,fH,S,lds,fh,h,V,matctx);
789: matctx->compM1 = PETSC_FALSE;
790: break;
791: case PEP_REFINE_SCHEME_SCHUR:
792: KSPGetOperators(ksp,&M,NULL);
793: MatShellGetContext(M,&ctx);
794: NRefSysSetup_shell(pep,k,fH,S,lds,fh,h,ctx);
795: ctx->compM1 = PETSC_FALSE;
796: break;
797: }
798: } else {
799: if (i%matctx->subc->n==0 && (idx=i+matctx->subc->color)<k) {
800: for (j=0;j<pep->nmat;j++) fh[j] = fH[j*k+idx+idx*lda];
801: h = H[idx+idx*ldh];
802: matctx->idx = idx;
803: switch (pep->scheme) {
804: case PEP_REFINE_SCHEME_EXPLICIT:
805: NRefSysSetup_explicit(pep,k,ksp,fH,S,lds,fh,h,matctx->V,matctx,matctx->W);
806: matctx->compM1 = PETSC_FALSE;
807: break;
808: case PEP_REFINE_SCHEME_MBE:
809: NRefSysSetup_mbe(pep,k,ksp,fH,S,lds,fh,h,matctx->V,matctx);
810: matctx->compM1 = PETSC_FALSE;
811: break;
812: case PEP_REFINE_SCHEME_SCHUR:
813: break;
814: }
815: } else idx = matctx->idx;
816: VecScatterBegin(matctx->scatter_id[i%matctx->subc->n],Rv,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
817: VecScatterEnd(matctx->scatter_id[i%matctx->subc->n],Rv,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
818: VecGetArrayRead(matctx->tg,&array);
819: VecPlaceArray(matctx->t,array);
820: VecCopy(matctx->t,matctx->Rv);
821: VecResetArray(matctx->t);
822: VecRestoreArrayRead(matctx->tg,&array);
823: R = matctx->Rv;
824: Vi = matctx->Vi;
825: }
826: if (idx==i && idx<k) {
827: switch (pep->scheme) {
828: case PEP_REFINE_SCHEME_EXPLICIT:
829: NRefSysSolve_explicit(k,ksp,R,Rh,Vi,dHi,matctx);
830: break;
831: case PEP_REFINE_SCHEME_MBE:
832: NRefSysSolve_mbe(k,k,matctx->W,matctx->w,matctx->Wt,matctx->wt,matctx->d,matctx->dt,ksp,matctx->M2,matctx->M3 ,matctx->M4,PETSC_FALSE,R,Rh,Vi,dHi,matctx->t);
833: break;
834: case PEP_REFINE_SCHEME_SCHUR:
835: NRefSysSolve_shell(ksp,pep->nmat,R,Rh,k,Vi,dHi);
836: break;
837: }
838: }
839: if (matctx && matctx->subc) {
840: VecGetLocalSize(Vi,&m);
841: VecGetArrayRead(Vi,&array);
842: VecGetArray(matctx->tg,&array2);
843: PetscMemcpy(array2,array,m*sizeof(PetscScalar));
844: VecRestoreArray(matctx->tg,&array2);
845: VecRestoreArrayRead(Vi,&array);
846: VecScatterBegin(matctx->scatter_id[i%matctx->subc->n],matctx->tg,dVi,INSERT_VALUES,SCATTER_REVERSE);
847: VecScatterEnd(matctx->scatter_id[i%matctx->subc->n],matctx->tg,dVi,INSERT_VALUES,SCATTER_REVERSE);
848: switch (pep->scheme) {
849: case PEP_REFINE_SCHEME_EXPLICIT:
850: MatGetOwnershipRange(matctx->E[0],&n0,&m0);
851: VecGetArrayRead(matctx->ttN,&array);
852: VecPlaceArray(matctx->tp,array+m0-n0);
853: VecScatterBegin(matctx->scatterp_id[i%matctx->subc->n],matctx->tp,matctx->tpg,INSERT_VALUES,SCATTER_FORWARD);
854: VecScatterEnd(matctx->scatterp_id[i%matctx->subc->n],matctx->tp,matctx->tpg,INSERT_VALUES,SCATTER_FORWARD);
855: VecResetArray(matctx->tp);
856: VecRestoreArrayRead(matctx->ttN,&array);
857: VecGetArrayRead(matctx->tpg,&array);
858: for (j=0;j<k;j++) dHi[j] = array[j];
859: VecRestoreArrayRead(matctx->tpg,&array);
860: break;
861: case PEP_REFINE_SCHEME_MBE:
862: root = 0;
863: for (j=0;j<i%matctx->subc->n;j++) root += matctx->subc->subsize[j];
864: PetscMPIIntCast(k,&len);
865: MPI_Bcast(dHi,len,MPIU_SCALAR,root,matctx->subc->dupparent);
866: break;
867: case PEP_REFINE_SCHEME_SCHUR:
868: break;
869: }
870: }
871: return(0);
872: }
876: static PetscErrorCode PEPNRefForwardSubstitution(PEP pep,PetscInt k,PetscScalar *S,PetscInt lds,PetscScalar *H,PetscInt ldh,PetscScalar *fH,BV dV,PetscScalar *dVS,PetscInt *rds,PetscScalar *dH,PetscInt lddh,KSP ksp,MatExplicitCtx *matctx)
877: {
879: PetscInt i,nmat=pep->nmat,lda=nmat*k;
880: PetscScalar *fh,*Rh,*DfH;
881: PetscReal norm;
882: BV W;
883: Vec Rv,t,dvi;
884: FSubctx *ctx;
885: Mat M,*At;
886: PetscBool flg,lindep;
889: PetscMalloc2(nmat*k*k,&DfH,k,&Rh);
890: *rds = 0;
891: BVCreateVec(pep->V,&Rv);
892: switch (pep->scheme) {
893: case PEP_REFINE_SCHEME_EXPLICIT:
894: BVCreateVec(pep->V,&t);
895: BVDuplicateResize(pep->V,PetscMax(k,nmat),&W);
896: PetscMalloc1(nmat,&fh);
897: break;
898: case PEP_REFINE_SCHEME_MBE:
899: if (matctx->subc) {
900: BVCreateVec(pep->V,&t);
901: BVDuplicateResize(pep->V,PetscMax(k,nmat),&W);
902: } else {
903: W = matctx->W;
904: PetscObjectReference((PetscObject)W);
905: t = matctx->t;
906: PetscObjectReference((PetscObject)t);
907: }
908: BVScale(matctx->W,0.0);
909: BVScale(matctx->Wt,0.0);
910: BVScale(matctx->M2,0.0);
911: BVScale(matctx->M3,0.0);
912: PetscMalloc1(nmat,&fh);
913: break;
914: case PEP_REFINE_SCHEME_SCHUR:
915: KSPGetOperators(ksp,&M,NULL);
916: MatShellGetContext(M,&ctx);
917: BVCreateVec(pep->V,&t);
918: BVDuplicateResize(pep->V,PetscMax(k,nmat),&W);
919: fh = ctx->fih;
920: break;
921: }
922: PetscMemzero(dVS,2*k*k*sizeof(PetscScalar));
923: PetscMemzero(DfH,lda*k*sizeof(PetscScalar));
924: STGetTransform(pep->st,&flg);
925: if (flg) {
926: PetscMalloc1(pep->nmat,&At);
927: for (i=0;i<pep->nmat;i++) {
928: STGetTOperators(pep->st,i,&At[i]);
929: }
930: } else At = pep->A;
932: /* Main loop for computing the ith columns of dX and dS */
933: for (i=0;i<k;i++) {
934: /* Compute and update i-th column of the right hand side */
935: PetscMemzero(Rh,k*sizeof(PetscScalar));
936: NRefRightSide(nmat,pep->pbc,At,k,pep->V,S,lds,i,H,ldh,fH,DfH,dH,dV,dVS,*rds,Rv,Rh,W,t);
938: /* Update and solve system */
939: BVGetColumn(dV,i,&dvi);
940: NRefSysIter(i,pep,k,ksp,fH,S,lds,fh,H,ldh,Rv,Rh,pep->V,dvi,dH+i*k,matctx,W);
941: /* Orthogonalize computed solution */
942: BVOrthogonalizeVec(pep->V,dvi,dVS+i*2*k,&norm,&lindep);
943: BVRestoreColumn(dV,i,&dvi);
944: if (!lindep) {
945: BVOrthogonalizeColumn(dV,i,dVS+k+i*2*k,&norm,&lindep);
946: if (!lindep) {
947: dVS[k+i+i*2*k] = norm;
948: BVScaleColumn(dV,i,1.0/norm);
949: (*rds)++;
950: }
951: }
952: }
953: BVSetActiveColumns(dV,0,*rds);
954: VecDestroy(&t);
955: VecDestroy(&Rv);
956: BVDestroy(&W);
957: if (flg) {
958: PetscFree(At);
959: }
960: PetscFree2(DfH,Rh);
961: if (pep->scheme!=PEP_REFINE_SCHEME_SCHUR) { PetscFree(fh); }
962: return(0);
963: }
967: static PetscErrorCode NRefOrthogStep(PEP pep,PetscInt k,PetscScalar *H,PetscInt ldh,PetscScalar *fH,PetscScalar *S,PetscInt lds,PetscInt *prs)
968: {
969: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
971: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
972: #else
974: PetscInt i,j,nmat=pep->nmat,deg=nmat-1,lda=nmat*k,rs=*prs,ldg;
975: PetscScalar *T,*G,*tau,*array,sone=1.0,zero=0.0,*work;
976: PetscBLASInt rs_,lds_,k_,ldh_,info,ldg_,lda_;
977: Mat M0;
980: PetscMalloc4(rs*k,&T,k,&tau,k,&work,deg*k*k,&G);
981: PetscBLASIntCast(lds,&lds_);
982: PetscBLASIntCast(lda,&lda_);
983: PetscBLASIntCast(k,&k_);
984: if (rs>k) { /* Truncate S to have k columns*/
985: for (j=0;j<k;j++) {
986: PetscMemcpy(T+j*rs,S+j*lds,rs*sizeof(PetscScalar));
987: }
988: PetscBLASIntCast(rs,&rs_);
989: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rs_,&k_,T,&rs_,tau,work,&k_,&info));
990: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
991: /* Copy triangular matrix in S */
992: PetscMemzero(S,lds*k*sizeof(PetscScalar));
993: for (j=0;j<k;j++) for (i=0;i<=j;i++) S[j*lds+i] = T[j*rs+i];
994: PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&rs_,&k_,&k_,T,&rs_,tau,work,&k_,&info));
995: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
996: MatCreateSeqDense(PETSC_COMM_SELF,rs,k,NULL,&M0);
997: MatDenseGetArray(M0,&array);
998: for (j=0;j<k;j++) {
999: PetscMemcpy(array+j*rs,T+j*rs,rs*sizeof(PetscScalar));
1000: }
1001: MatDenseRestoreArray(M0,&array);
1002: BVSetActiveColumns(pep->V,0,rs);
1003: BVMultInPlace(pep->V,M0,0,k);
1004: BVSetActiveColumns(pep->V,0,k);
1005: MatDestroy(&M0);
1006: *prs = rs = k;
1007: }
1008: /* Form auxiliary matrix for the orthogonalization step */
1009: ldg = deg*k;
1010: PEPEvaluateBasisforMatrix(pep,nmat,k,H,ldh,fH);
1011: PetscBLASIntCast(ldg,&ldg_);
1012: PetscBLASIntCast(ldh,&ldh_);
1013: for (j=0;j<deg;j++) {
1014: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&sone,S,&lds_,fH+j*k,&lda_,&zero,G+j*k,&ldg_));
1015: }
1016: /* Orthogonalize and update S */
1017: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&ldg_,&k_,G,&ldg_,tau,work,&k_,&info));
1018: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
1019: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&k_,&k_,&sone,G,&ldg_,S,&lds_));
1021: /* Update H */
1022: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&k_,&k_,&sone,G,&ldg_,H,&ldh_));
1023: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("R","U","N","N",&k_,&k_,&sone,G,&ldg_,H,&ldh_));
1024: PetscFree4(T,tau,work,G);
1025: return(0);
1026: #endif
1027: }
1031: static PetscErrorCode PEPNRefUpdateInvPair(PEP pep,PetscInt k,PetscScalar *H,PetscInt ldh,PetscScalar *fH,PetscScalar *dH,PetscScalar *S,PetscInt lds,BV dV,PetscScalar *dVS,PetscInt rds)
1032: {
1033: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
1035: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
1036: #else
1038: PetscInt i,j,nmat=pep->nmat,lda=nmat*k;
1039: PetscScalar *tau,*array,*work;
1040: PetscBLASInt lds_,k_,lda_,ldh_,kdrs_,info,k2_;
1041: Mat M0;
1044: PetscMalloc2(k,&tau,k,&work);
1045: PetscBLASIntCast(lds,&lds_);
1046: PetscBLASIntCast(lda,&lda_);
1047: PetscBLASIntCast(ldh,&ldh_);
1048: PetscBLASIntCast(k,&k_);
1049: PetscBLASIntCast(2*k,&k2_);
1050: PetscBLASIntCast((k+rds),&kdrs_);
1051: /* Update H */
1052: for (j=0;j<k;j++) {
1053: for (i=0;i<k;i++) H[i+j*ldh] -= dH[i+j*k];
1054: }
1055: /* Update V */
1056: for (j=0;j<k;j++) {
1057: for (i=0;i<k;i++) dVS[i+j*2*k] = -dVS[i+j*2*k]+S[i+j*lds];
1058: for (i=k;i<2*k;i++) dVS[i+j*2*k] = -dVS[i+j*2*k];
1059: }
1060: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&kdrs_,&k_,dVS,&k2_,tau,work,&k_,&info));
1061: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
1062: /* Copy triangular matrix in S */
1063: for (j=0;j<k;j++) {
1064: for (i=0;i<=j;i++) S[i+j*lds] = dVS[i+j*2*k];
1065: for (i=j+1;i<k;i++) S[i+j*lds] = 0.0;
1066: }
1067: PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&k2_,&k_,&k_,dVS,&k2_,tau,work,&k_,&info));
1068: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
1069: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M0);
1070: MatDenseGetArray(M0,&array);
1071: for (j=0;j<k;j++) {
1072: PetscMemcpy(array+j*k,dVS+j*2*k,k*sizeof(PetscScalar));
1073: }
1074: MatDenseRestoreArray(M0,&array);
1075: BVMultInPlace(pep->V,M0,0,k);
1076: if (rds) {
1077: MatDenseGetArray(M0,&array);
1078: for (j=0;j<k;j++) {
1079: PetscMemcpy(array+j*k,dVS+k+j*2*k,rds*sizeof(PetscScalar));
1080: }
1081: MatDenseRestoreArray(M0,&array);
1082: BVMultInPlace(dV,M0,0,k);
1083: BVMult(pep->V,1.0,1.0,dV,NULL);
1084: }
1085: MatDestroy(&M0);
1086: NRefOrthogStep(pep,k,H,ldh,fH,S,lds,&k);
1087: PetscFree2(tau,work);
1088: return(0);
1089: #endif
1090: }
1094: static PetscErrorCode PEPNRefSetUp(PEP pep,PetscInt k,PetscScalar *H,PetscInt ldh,MatExplicitCtx *matctx,PetscBool ini)
1095: {
1096: PetscErrorCode ierr;
1097: FSubctx *ctx;
1098: PetscScalar t,*coef;
1099: const PetscScalar *array;
1100: MatStructure str;
1101: PetscInt j,nmat=pep->nmat,n0,m0,n1,m1,n0_,m0_,n1_,m1_,N0,N1,p,*idx1,*idx2,count,si,i,l0;
1102: MPI_Comm comm;
1103: PetscMPIInt np;
1104: const PetscInt *rgs0,*rgs1;
1105: Mat B,C,*E,*A,*At;
1106: IS is1,is2;
1107: Vec v;
1108: PetscBool flg;
1109: Mat M,P;
1112: PetscMalloc1(nmat,&coef);
1113: STGetTransform(pep->st,&flg);
1114: if (flg) {
1115: PetscMalloc1(pep->nmat,&At);
1116: for (i=0;i<pep->nmat;i++) {
1117: STGetTOperators(pep->st,i,&At[i]);
1118: }
1119: } else At = pep->A;
1120: switch (pep->scheme) {
1121: case PEP_REFINE_SCHEME_EXPLICIT:
1122: if (ini) {
1123: if (matctx->subc) {
1124: A = matctx->A;
1125: comm = PetscSubcommChild(matctx->subc);
1126: } else {
1127: A = At;
1128: PetscObjectGetComm((PetscObject)pep,&comm);
1129: }
1130: E = matctx->E;
1131: STGetMatStructure(pep->st,&str);
1132: MatDuplicate(A[0],MAT_COPY_VALUES,&E[0]);
1133: j = (matctx->subc)?matctx->subc->color:0;
1134: PEPEvaluateBasis(pep,H[j+j*ldh],0,coef,NULL);
1135: for (j=1;j<nmat;j++) {
1136: MatAXPY(E[0],coef[j],A[j],str);
1137: }
1138: MatCreateDense(comm,PETSC_DECIDE,PETSC_DECIDE,k,k,NULL,&E[1]);
1139: MatAssemblyBegin(E[1],MAT_FINAL_ASSEMBLY);
1140: MatAssemblyEnd(E[1],MAT_FINAL_ASSEMBLY);
1141: MatGetOwnershipRange(E[0],&n0,&m0);
1142: MatGetOwnershipRange(E[1],&n1,&m1);
1143: MatGetOwnershipRangeColumn(E[0],&n0_,&m0_);
1144: MatGetOwnershipRangeColumn(E[1],&n1_,&m1_);
1145: /* T12 and T21 are computed from V and V*, so,
1146: they must have the same column and row ranges */
1147: if (m0_-n0_ != m0-n0) SETERRQ(PETSC_COMM_SELF,1,"Inconsistent dimensions");
1148: MatCreateDense(comm,m0-n0,m1_-n1_,PETSC_DECIDE,PETSC_DECIDE,NULL,&B);
1149: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1150: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1151: MatCreateDense(comm,m1-n1,m0_-n0_,PETSC_DECIDE,PETSC_DECIDE,NULL,&C);
1152: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1153: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1154: SlepcMatTile(1.0,E[0],1.0,B,1.0,C,1.0,E[1],&M);
1155: MatDestroy(&B);
1156: MatDestroy(&C);
1157: matctx->compM1 = PETSC_TRUE;
1158: MatGetSize(E[0],NULL,&N0);
1159: MatGetSize(E[1],NULL,&N1);
1160: MPI_Comm_size(PetscObjectComm((PetscObject)M),&np);
1161: MatGetOwnershipRanges(E[0],&rgs0);
1162: MatGetOwnershipRanges(E[1],&rgs1);
1163: PetscMalloc4(PetscMax(k,N1),&matctx->idxp,N0,&matctx->idxg,N0,&matctx->map0,N1,&matctx->map1);
1164: /* Create column (and row) mapping */
1165: for (p=0;p<np;p++) {
1166: for (j=rgs0[p];j<rgs0[p+1];j++) matctx->map0[j] = j+rgs1[p];
1167: for (j=rgs1[p];j<rgs1[p+1];j++) matctx->map1[j] = j+rgs0[p+1];
1168: }
1169: MatCreateVecs(M,NULL,&matctx->tN);
1170: MatCreateVecs(matctx->E[1],NULL,&matctx->t1);
1171: VecDuplicate(matctx->tN,&matctx->ttN);
1172: if (matctx->subc) {
1173: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1174: count = np*k;
1175: PetscMalloc2(count,&idx1,count,&idx2);
1176: VecCreateMPI(PetscObjectComm((PetscObject)pep),m1-n1,PETSC_DECIDE,&matctx->tp);
1177: VecGetOwnershipRange(matctx->tp,&l0,NULL);
1178: VecCreateMPI(PetscObjectComm((PetscObject)pep),k,PETSC_DECIDE,&matctx->tpg);
1179: for (si=0;si<matctx->subc->n;si++) {
1180: if (matctx->subc->color==si) {
1181: j=0;
1182: if (matctx->subc->color==si) {
1183: for (p=0;p<np;p++) {
1184: for (i=n1;i<m1;i++) {
1185: idx1[j] = l0+i-n1;
1186: idx2[j++] =p*k+i;
1187: }
1188: }
1189: }
1190: count = np*(m1-n1);
1191: } else count =0;
1192: ISCreateGeneral(PetscObjectComm((PetscObject)pep),count,idx1,PETSC_COPY_VALUES,&is1);
1193: ISCreateGeneral(PetscObjectComm((PetscObject)pep),count,idx2,PETSC_COPY_VALUES,&is2);
1194: VecScatterCreate(matctx->tp,is1,matctx->tpg,is2,&matctx->scatterp_id[si]);
1195: ISDestroy(&is1);
1196: ISDestroy(&is2);
1197: }
1198: PetscFree2(idx1,idx2);
1199: } else {
1200: VecScatterCreateToAll(matctx->t1,&matctx->scatterctx,&matctx->vseq);
1201: }
1202: P = M;
1203: } else {
1204: if (matctx->subc) {
1205: /* Scatter vectors pep->V */
1206: for (i=0;i<k;i++) {
1207: BVGetColumn(pep->V,i,&v);
1208: VecScatterBegin(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1209: VecScatterEnd(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1210: BVRestoreColumn(pep->V,i,&v);
1211: VecGetArrayRead(matctx->tg,&array);
1212: VecPlaceArray(matctx->t,(const PetscScalar*)array);
1213: BVInsertVec(matctx->V,i,matctx->t);
1214: VecResetArray(matctx->t);
1215: VecRestoreArrayRead(matctx->tg,&array);
1216: }
1217: }
1218: }
1219: break;
1220: case PEP_REFINE_SCHEME_MBE:
1221: if (ini) {
1222: if (matctx->subc) {
1223: A = matctx->A;
1224: comm = PetscSubcommChild(matctx->subc);
1225: } else {
1226: matctx->V = pep->V;
1227: A = At;
1228: PetscObjectGetComm((PetscObject)pep,&comm);
1229: MatCreateVecs(pep->A[0],&matctx->t,NULL);
1230: }
1231: STGetMatStructure(pep->st,&str);
1232: MatDuplicate(A[0],MAT_COPY_VALUES,&matctx->M1);
1233: j = (matctx->subc)?matctx->subc->color:0;
1234: PEPEvaluateBasis(pep,H[j+j*ldh],0,coef,NULL);
1235: for (j=1;j<nmat;j++) {
1236: MatAXPY(matctx->M1,coef[j],A[j],str);
1237: }
1238: BVDuplicateResize(matctx->V,PetscMax(k,pep->nmat),&matctx->W);
1239: BVDuplicateResize(matctx->V,k,&matctx->M2);
1240: BVDuplicate(matctx->M2,&matctx->M3);
1241: BVDuplicate(matctx->M2,&matctx->Wt);
1242: PetscMalloc5(k*k,&matctx->M4,k*k,&matctx->w,k*k,&matctx->wt,k,&matctx->d,k,&matctx->dt);
1243: matctx->compM1 = PETSC_TRUE;
1244: M = matctx->M1;
1245: P = M;
1246: }
1247: break;
1248: case PEP_REFINE_SCHEME_SCHUR:
1249: if (ini) {
1250: PetscObjectGetComm((PetscObject)pep,&comm);
1251: MatGetSize(At[0],&m0,&n0);
1252: PetscMalloc1(1,&ctx);
1253: STGetMatStructure(pep->st,&str);
1254: /* Create a shell matrix to solve the linear system */
1255: ctx->A = At;
1256: ctx->V = pep->V;
1257: ctx->k = k; ctx->nmat = nmat;
1258: PetscMalloc4(k*k,&ctx->M4,k,&ctx->pM4,2*k*k,&ctx->work,nmat,&ctx->fih);
1259: PetscMemzero(ctx->M4,k*k*sizeof(PetscScalar));
1260: MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,m0,n0,ctx,&M);
1261: MatShellSetOperation(M,MATOP_MULT,(void(*)(void))MatFSMult);
1262: BVDuplicateResize(ctx->V,PetscMax(k,pep->nmat),&ctx->W);
1263: BVDuplicateResize(ctx->V,k,&ctx->M2);
1264: BVDuplicate(ctx->M2,&ctx->M3);
1265: BVCreateVec(pep->V,&ctx->t);
1266: MatDuplicate(At[0],MAT_COPY_VALUES,&ctx->M1);
1267: PEPEvaluateBasis(pep,H[0],0,coef,NULL);
1268: for (j=1;j<nmat;j++) {
1269: MatAXPY(ctx->M1,coef[j],At[j],str);
1270: }
1271: MatDuplicate(At[0],MAT_COPY_VALUES,&P);
1272: /* Compute a precond matrix for the system */
1273: t = H[0];
1274: PEPEvaluateBasis(pep,t,0,coef,NULL);
1275: for (j=1;j<nmat;j++) {
1276: MatAXPY(P,coef[j],At[j],str);
1277: }
1278: ctx->compM1 = PETSC_TRUE;
1279: }
1280: break;
1281: }
1282: if (ini) {
1283: PEPRefineGetKSP(pep,&pep->refineksp);
1284: KSPSetErrorIfNotConverged(pep->refineksp,PETSC_TRUE);
1285: KSPSetOperators(pep->refineksp,M,P);
1286: KSPSetFromOptions(pep->refineksp);
1287: }
1289: if (!ini && matctx && matctx->subc) {
1290: /* Scatter vectors pep->V */
1291: for (i=0;i<k;i++) {
1292: BVGetColumn(pep->V,i,&v);
1293: VecScatterBegin(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1294: VecScatterEnd(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1295: BVRestoreColumn(pep->V,i,&v);
1296: VecGetArrayRead(matctx->tg,&array);
1297: VecPlaceArray(matctx->t,(const PetscScalar*)array);
1298: BVInsertVec(matctx->V,i,matctx->t);
1299: VecResetArray(matctx->t);
1300: VecRestoreArrayRead(matctx->tg,&array);
1301: }
1302: }
1303: PetscFree(coef);
1304: if (flg) {
1305: PetscFree(At);
1306: }
1307: return(0);
1308: }
1312: static PetscErrorCode NRefSubcommSetup(PEP pep,PetscInt k,MatExplicitCtx *matctx,PetscInt nsubc)
1313: {
1314: PetscErrorCode ierr;
1315: PetscInt i,si,j,m0,n0,nloc0,nloc_sub,*idx1,*idx2;
1316: IS is1,is2;
1317: BVType type;
1318: Vec v;
1319: const PetscScalar *array;
1320: Mat *A;
1321: PetscBool flg;
1324: STGetTransform(pep->st,&flg);
1325: if (flg) {
1326: PetscMalloc1(pep->nmat,&A);
1327: for (i=0;i<pep->nmat;i++) {
1328: STGetTOperators(pep->st,i,&A[i]);
1329: }
1330: } else A = pep->A;
1332: /* Duplicate pep matrices */
1333: PetscMalloc3(pep->nmat,&matctx->A,nsubc,&matctx->scatter_id,nsubc,&matctx->scatterp_id);
1334: for (i=0;i<pep->nmat;i++) {
1335: MatCreateRedundantMatrix(A[i],0,PetscSubcommChild(matctx->subc),MAT_INITIAL_MATRIX,&matctx->A[i]);
1336: }
1338: /* Create Scatter */
1339: MatCreateVecs(matctx->A[0],&matctx->t,NULL);
1340: MatGetLocalSize(matctx->A[0],&nloc_sub,NULL);
1341: VecCreateMPI(PetscSubcommContiguousParent(matctx->subc),nloc_sub,PETSC_DECIDE,&matctx->tg);
1342: BVGetColumn(pep->V,0,&v);
1343: VecGetOwnershipRange(v,&n0,&m0);
1344: nloc0 = m0-n0;
1345: PetscMalloc2(matctx->subc->n*nloc0,&idx1,matctx->subc->n*nloc0,&idx2);
1346: j = 0;
1347: for (si=0;si<matctx->subc->n;si++) {
1348: for (i=n0;i<m0;i++) {
1349: idx1[j] = i;
1350: idx2[j++] = i+pep->n*si;
1351: }
1352: }
1353: ISCreateGeneral(PetscObjectComm((PetscObject)pep),matctx->subc->n*nloc0,idx1,PETSC_COPY_VALUES,&is1);
1354: ISCreateGeneral(PetscObjectComm((PetscObject)pep),matctx->subc->n*nloc0,idx2,PETSC_COPY_VALUES,&is2);
1355: VecScatterCreate(v,is1,matctx->tg,is2,&matctx->scatter_sub);
1356: ISDestroy(&is1);
1357: ISDestroy(&is2);
1358: for (si=0;si<matctx->subc->n;si++) {
1359: j=0;
1360: for (i=n0;i<m0;i++) {
1361: idx1[j] = i;
1362: idx2[j++] = i+pep->n*si;
1363: }
1364: ISCreateGeneral(PetscObjectComm((PetscObject)pep),nloc0,idx1,PETSC_COPY_VALUES,&is1);
1365: ISCreateGeneral(PetscObjectComm((PetscObject)pep),nloc0,idx2,PETSC_COPY_VALUES,&is2);
1366: VecScatterCreate(v,is1,matctx->tg,is2,&matctx->scatter_id[si]);
1367: ISDestroy(&is1);
1368: ISDestroy(&is2);
1369: }
1370: BVRestoreColumn(pep->V,0,&v);
1371: PetscFree2(idx1,idx2);
1373: /* Duplicate pep->V vecs */
1374: BVGetType(pep->V,&type);
1375: BVCreate(PetscSubcommChild(matctx->subc),&matctx->V);
1376: BVSetType(matctx->V,type);
1377: BVSetSizesFromVec(matctx->V,matctx->t,k);
1378: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT) {
1379: BVDuplicateResize(matctx->V,PetscMax(k,pep->nmat),&matctx->W);
1380: }
1381: for (i=0;i<k;i++) {
1382: BVGetColumn(pep->V,i,&v);
1383: VecScatterBegin(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1384: VecScatterEnd(matctx->scatter_sub,v,matctx->tg,INSERT_VALUES,SCATTER_FORWARD);
1385: BVRestoreColumn(pep->V,i,&v);
1386: VecGetArrayRead(matctx->tg,&array);
1387: VecPlaceArray(matctx->t,(const PetscScalar*)array);
1388: BVInsertVec(matctx->V,i,matctx->t);
1389: VecResetArray(matctx->t);
1390: VecRestoreArrayRead(matctx->tg,&array);
1391: }
1393: VecDuplicate(matctx->t,&matctx->Rv);
1394: VecDuplicate(matctx->t,&matctx->Vi);
1395: if (flg) {
1396: PetscFree(A);
1397: }
1398: return(0);
1399: }
1403: static PetscErrorCode NRefSubcommDestroy(PEP pep,MatExplicitCtx *matctx)
1404: {
1406: PetscInt i;
1409: VecScatterDestroy(&matctx->scatter_sub);
1410: for (i=0;i<matctx->subc->n;i++) {
1411: VecScatterDestroy(&matctx->scatter_id[i]);
1412: }
1413: for (i=0;i<pep->nmat;i++) {
1414: MatDestroy(&matctx->A[i]);
1415: }
1416: if (pep->scheme==PEP_REFINE_SCHEME_EXPLICIT) {
1417: for (i=0;i<matctx->subc->n;i++) {
1418: VecScatterDestroy(&matctx->scatterp_id[i]);
1419: }
1420: VecDestroy(&matctx->tp);
1421: VecDestroy(&matctx->tpg);
1422: BVDestroy(&matctx->W);
1423: }
1424: PetscFree3(matctx->A,matctx->scatter_id,matctx->scatterp_id);
1425: BVDestroy(&matctx->V);
1426: VecDestroy(&matctx->t);
1427: VecDestroy(&matctx->tg);
1428: VecDestroy(&matctx->Rv);
1429: VecDestroy(&matctx->Vi);
1430: return(0);
1431: }
1435: PetscErrorCode PEPNewtonRefinement_TOAR(PEP pep,PetscScalar sigma,PetscInt *maxits,PetscReal *tol,PetscInt k,PetscScalar *S,PetscInt lds,PetscInt *prs)
1436: {
1437: #if defined(PETSC_MISSING_LAPACK_GETRF) || defined(PETSC_MISSING_LAPACK_GETRI)
1439: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF/GETRI - Lapack routine is unavailable");
1440: #else
1442: PetscScalar *H,*work,*dH,*fH,*dVS;
1443: PetscInt ldh,i,j,its=1,nmat=pep->nmat,nsubc=pep->npart,rds;
1444: PetscBLASInt k_,ld_,*p,info;
1445: BV dV;
1446: PetscBool sinvert,flg;
1447: MatExplicitCtx *matctx=NULL;
1448: Vec v;
1449: Mat M,P;
1450: FSubctx *ctx;
1453: PetscLogEventBegin(PEP_Refine,pep,0,0,0);
1454: if (k > pep->n) SETERRQ1(PetscObjectComm((PetscObject)pep),1,"Multiple Refinement available only for invariant pairs of dimension smaller than n=%D",pep->n);
1455: /* the input tolerance is not being taken into account (by the moment) */
1456: its = *maxits;
1457: PetscMalloc3(k*k,&dH,nmat*k*k,&fH,k,&work);
1458: DSGetLeadingDimension(pep->ds,&ldh);
1459: DSGetArray(pep->ds,DS_MAT_A,&H);
1460: DSRestoreArray(pep->ds,DS_MAT_A,&H);
1461: PetscMalloc1(2*k*k,&dVS);
1462: STGetTransform(pep->st,&flg);
1463: if (!flg && pep->st && pep->ops->backtransform) { /* BackTransform */
1464: PetscBLASIntCast(k,&k_);
1465: PetscBLASIntCast(ldh,&ld_);
1466: PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert);
1467: if (sinvert) {
1468: DSGetArray(pep->ds,DS_MAT_A,&H);
1469: PetscMalloc1(k,&p);
1470: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,H,&ld_,p,&info));
1471: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,H,&ld_,p,work,&k_,&info));
1472: DSRestoreArray(pep->ds,DS_MAT_A,&H);
1473: pep->ops->backtransform = NULL;
1474: }
1475: if (sigma!=0.0) {
1476: DSGetArray(pep->ds,DS_MAT_A,&H);
1477: for (i=0;i<k;i++) H[i+ldh*i] += sigma;
1478: DSRestoreArray(pep->ds,DS_MAT_A,&H);
1479: pep->ops->backtransform = NULL;
1480: }
1481: }
1482: if ((pep->scale==PEP_SCALE_BOTH || pep->scale==PEP_SCALE_SCALAR) && pep->sfactor!=1.0) {
1483: DSGetArray(pep->ds,DS_MAT_A,&H);
1484: for (j=0;j<k;j++) {
1485: for (i=0;i<k;i++) H[i+j*ldh] *= pep->sfactor;
1486: }
1487: DSRestoreArray(pep->ds,DS_MAT_A,&H);
1488: if (!flg) {
1489: /* Restore original values */
1490: for (i=0;i<pep->nmat;i++){
1491: pep->pbc[pep->nmat+i] *= pep->sfactor;
1492: pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
1493: }
1494: }
1495: }
1496: if ((pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) && pep->Dr) {
1497: for (i=0;i<k;i++) {
1498: BVGetColumn(pep->V,i,&v);
1499: VecPointwiseMult(v,v,pep->Dr);
1500: BVRestoreColumn(pep->V,i,&v);
1501: }
1502: }
1503: DSGetArray(pep->ds,DS_MAT_A,&H);
1505: NRefOrthogStep(pep,k,H,ldh,fH,S,lds,prs);
1506: /* check if H is in Schur form */
1507: for (i=0;i<k-1;i++) {
1508: if (H[i+1+i*ldh]!=0.0) {
1509: #if !defined(PETSC_USES_COMPLEX)
1510: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Iterative Refinement requires the complex Schur form of the projected matrix");
1511: #else
1512: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Iterative Refinement requires an upper triangular projected matrix");
1513: #endif
1514: }
1515: }
1516: if (nsubc>k) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Amount of subcommunicators should not be larger than the invariant pair dimension");
1517: BVSetActiveColumns(pep->V,0,k);
1518: BVDuplicateResize(pep->V,k,&dV);
1519: PetscLogObjectParent((PetscObject)pep,(PetscObject)dV);
1520: if (pep->scheme!=PEP_REFINE_SCHEME_SCHUR) {
1521: PetscMalloc1(1,&matctx);
1522: if (nsubc>1) { /* spliting in subcommunicators */
1523: matctx->subc = pep->refinesubc;
1524: NRefSubcommSetup(pep,k,matctx,nsubc);
1525: } else matctx->subc=NULL;
1526: }
1528: /* Loop performing iterative refinements */
1529: for (i=0;i<its;i++) {
1530: /* Pre-compute the polynomial basis evaluated in H */
1531: PEPEvaluateBasisforMatrix(pep,nmat,k,H,ldh,fH);
1532: PEPNRefSetUp(pep,k,H,ldh,matctx,PetscNot(i));
1533: /* Solve the linear system */
1534: PEPNRefForwardSubstitution(pep,k,S,lds,H,ldh,fH,dV,dVS,&rds,dH,k,pep->refineksp,matctx);
1535: /* Update X (=V*S) and H, and orthogonalize [X;X*fH1;...;XfH(deg-1)] */
1536: PEPNRefUpdateInvPair(pep,k,H,ldh,fH,dH,S,lds,dV,dVS,rds);
1537: }
1538: DSRestoreArray(pep->ds,DS_MAT_A,&H);
1539: if (!flg && sinvert) {
1540: PetscFree(p);
1541: }
1542: PetscFree3(dH,fH,work);
1543: PetscFree(dVS);
1544: BVDestroy(&dV);
1545: switch (pep->scheme) {
1546: case PEP_REFINE_SCHEME_EXPLICIT:
1547: for (i=0;i<2;i++) {
1548: MatDestroy(&matctx->E[i]);
1549: }
1550: PetscFree4(matctx->idxp,matctx->idxg,matctx->map0,matctx->map1);
1551: VecDestroy(&matctx->tN);
1552: VecDestroy(&matctx->ttN);
1553: VecDestroy(&matctx->t1);
1554: if (nsubc>1) {
1555: NRefSubcommDestroy(pep,matctx);
1556: } else {
1557: VecDestroy(&matctx->vseq);
1558: VecScatterDestroy(&matctx->scatterctx);
1559: }
1560: PetscFree(matctx);
1561: KSPGetOperators(pep->refineksp,&M,NULL);
1562: MatDestroy(&M);
1563: break;
1564: case PEP_REFINE_SCHEME_MBE:
1565: BVDestroy(&matctx->W);
1566: BVDestroy(&matctx->Wt);
1567: BVDestroy(&matctx->M2);
1568: BVDestroy(&matctx->M3);
1569: MatDestroy(&matctx->M1);
1570: VecDestroy(&matctx->t);
1571: PetscFree5(matctx->M4,matctx->w,matctx->wt,matctx->d,matctx->dt);
1572: if (nsubc>1) {
1573: NRefSubcommDestroy(pep,matctx);
1574: }
1575: PetscFree(matctx);
1576: break;
1577: case PEP_REFINE_SCHEME_SCHUR:
1578: KSPGetOperators(pep->refineksp,&M,&P);
1579: MatShellGetContext(M,&ctx);
1580: PetscFree4(ctx->M4,ctx->pM4,ctx->work,ctx->fih);
1581: MatDestroy(&ctx->M1);
1582: BVDestroy(&ctx->M2);
1583: BVDestroy(&ctx->M3);
1584: BVDestroy(&ctx->W);
1585: VecDestroy(&ctx->t);
1586: PetscFree(ctx);
1587: MatDestroy(&M);
1588: MatDestroy(&P);
1589: break;
1590: }
1591: PetscLogEventEnd(PEP_Refine,pep,0,0,0);
1592: return(0);
1593: #endif
1594: }