Actual source code: blzpack.c
slepc-3.12.2 2020-01-13
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: This file implements a wrapper to the BLZPACK package
12: */
14: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
15: #include "blzpack.h"
17: const char* blzpack_error[33] = {
18: "",
19: "illegal data, LFLAG ",
20: "illegal data, dimension of (U), (V), (X) ",
21: "illegal data, leading dimension of (U), (V), (X) ",
22: "illegal data, leading dimension of (EIG) ",
23: "illegal data, number of required eigenpairs ",
24: "illegal data, Lanczos algorithm block size ",
25: "illegal data, maximum number of steps ",
26: "illegal data, number of starting vectors ",
27: "illegal data, number of eigenpairs provided ",
28: "illegal data, problem type flag ",
29: "illegal data, spectrum slicing flag ",
30: "illegal data, eigenvectors purification flag ",
31: "illegal data, level of output ",
32: "illegal data, output file unit ",
33: "illegal data, LCOMM (MPI or PVM) ",
34: "illegal data, dimension of ISTOR ",
35: "illegal data, convergence threshold ",
36: "illegal data, dimension of RSTOR ",
37: "illegal data on at least one PE ",
38: "ISTOR(3:14) must be equal on all PEs ",
39: "RSTOR(1:3) must be equal on all PEs ",
40: "not enough space in ISTOR to start eigensolution ",
41: "not enough space in RSTOR to start eigensolution ",
42: "illegal data, number of negative eigenvalues ",
43: "illegal data, entries of V ",
44: "illegal data, entries of X ",
45: "failure in computational subinterval ",
46: "file I/O error, blzpack.__.BQ ",
47: "file I/O error, blzpack.__.BX ",
48: "file I/O error, blzpack.__.Q ",
49: "file I/O error, blzpack.__.X ",
50: "parallel interface error "
51: };
53: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
54: {
56: PetscInt listor,lrstor,ncuv,k1,k2,k3,k4;
57: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
58: PetscBool issinv,istrivial,flg;
61: if (eps->ncv) {
62: if (eps->ncv < PetscMin(eps->nev+10,eps->nev*2)) SETERRQ(PetscObjectComm((PetscObject)eps),0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
63: } else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
64: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
65: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
67: if (!blz->block_size) blz->block_size = 3;
68: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
69: if (eps->which==EPS_ALL) {
70: if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Must define a computational interval when using EPS_ALL");
71: blz->slice = 1;
72: }
73: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&issinv);
74: if (blz->slice || eps->isgeneralized) {
75: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert ST is needed for generalized problems or spectrum slicing");
76: }
77: if (blz->slice) {
78: if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
79: if (eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The defined computational interval should have at least one of their sides bounded");
80: STSetDefaultShift(eps->st,eps->inta);
81: } else {
82: STSetDefaultShift(eps->st,eps->intb);
83: }
84: }
85: if (!eps->which) {
86: if (issinv) eps->which = EPS_TARGET_REAL;
87: else eps->which = EPS_SMALLEST_REAL;
88: }
89: if ((issinv && eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_MAGNITUDE && eps->which!=EPS_ALL) || (!issinv && eps->which!=EPS_SMALLEST_REAL)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
90: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
92: k1 = PetscMin(eps->n,180);
93: k2 = blz->block_size;
94: k4 = PetscMin(eps->ncv,eps->n);
95: k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;
97: listor = 123+k1*12;
98: PetscFree(blz->istor);
99: PetscMalloc1((17+listor),&blz->istor);
100: PetscLogObjectMemory((PetscObject)eps,(17+listor)*sizeof(PetscBLASInt));
101: PetscBLASIntCast(listor,&blz->istor[14]);
103: if (blz->slice) lrstor = eps->nloc*(k2*4+k1*2+k4)+k3;
104: else lrstor = eps->nloc*(k2*4+k1)+k3;
105: lrstor*=10;
106: PetscFree(blz->rstor);
107: PetscMalloc1((4+lrstor),&blz->rstor);
108: PetscLogObjectMemory((PetscObject)eps,(4+lrstor)*sizeof(PetscReal));
109: blz->rstor[3] = lrstor;
111: ncuv = PetscMax(3,blz->block_size);
112: PetscFree(blz->u);
113: PetscMalloc1(ncuv*eps->nloc,&blz->u);
114: PetscLogObjectMemory((PetscObject)eps,ncuv*eps->nloc*sizeof(PetscScalar));
115: PetscFree(blz->v);
116: PetscMalloc1(ncuv*eps->nloc,&blz->v);
117: PetscLogObjectMemory((PetscObject)eps,ncuv*eps->nloc*sizeof(PetscScalar));
119: PetscFree(blz->eig);
120: PetscMalloc1(2*eps->ncv,&blz->eig);
121: PetscLogObjectMemory((PetscObject)eps,2*eps->ncv*sizeof(PetscReal));
123: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
125: EPSAllocateSolution(eps,0);
126: EPS_SetInnerProduct(eps);
127: PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
128: if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
129: RGIsTrivial(eps->rg,&istrivial);
130: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
131: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
132: return(0);
133: }
135: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
136: {
138: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
139: PetscInt nn;
140: PetscBLASInt i,nneig,lflag,nvopu;
141: Vec x,y;
142: PetscScalar sigma,*pV;
143: Mat A,M;
144: KSP ksp;
145: PC pc;
148: STGetMatrix(eps->st,0,&A);
149: MatCreateVecsEmpty(A,&x,&y);
150: EPSGetStartVector(eps,0,NULL);
151: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
152: BVGetArray(eps->V,&pV);
154: if (eps->isgeneralized && !blz->slice) {
155: STGetShift(eps->st,&sigma); /* shift of origin */
156: blz->rstor[0] = sigma; /* lower limit of eigenvalue interval */
157: blz->rstor[1] = sigma; /* upper limit of eigenvalue interval */
158: } else {
159: sigma = 0.0;
160: blz->rstor[0] = eps->inta; /* lower limit of eigenvalue interval */
161: blz->rstor[1] = eps->intb; /* upper limit of eigenvalue interval */
162: }
163: nneig = 0; /* no. of eigs less than sigma */
165: PetscBLASIntCast(eps->nloc,&blz->istor[0]); /* no. of rows of U, V, X */
166: PetscBLASIntCast(eps->nloc,&blz->istor[1]); /* leading dim of U, V, X */
167: PetscBLASIntCast(eps->nev,&blz->istor[2]); /* required eigenpairs */
168: PetscBLASIntCast(eps->ncv,&blz->istor[3]); /* working eigenpairs */
169: blz->istor[4] = blz->block_size; /* number of vectors in a block */
170: blz->istor[5] = blz->nsteps; /* maximun number of steps per run */
171: blz->istor[6] = 1; /* number of starting vectors as input */
172: blz->istor[7] = 0; /* number of eigenpairs given as input */
173: blz->istor[8] = (blz->slice || eps->isgeneralized) ? 1 : 0; /* problem type */
174: blz->istor[9] = blz->slice; /* spectrum slicing */
175: blz->istor[10] = eps->isgeneralized ? 1 : 0; /* solutions refinement (purify) */
176: blz->istor[11] = 0; /* level of printing */
177: blz->istor[12] = 6; /* file unit for output */
178: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&blz->istor[13]);
180: blz->rstor[2] = eps->tol; /* threshold for convergence */
182: lflag = 0; /* reverse communication interface flag */
184: do {
185: BLZpack_(blz->istor,blz->rstor,&sigma,&nneig,blz->u,blz->v,&lflag,&nvopu,blz->eig,pV);
187: switch (lflag) {
188: case 1:
189: /* compute v = OP u */
190: for (i=0;i<nvopu;i++) {
191: VecPlaceArray(x,blz->u+i*eps->nloc);
192: VecPlaceArray(y,blz->v+i*eps->nloc);
193: if (blz->slice || eps->isgeneralized) {
194: STMatSolve(eps->st,x,y);
195: } else {
196: STApply(eps->st,x,y);
197: }
198: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
199: VecResetArray(x);
200: VecResetArray(y);
201: }
202: /* monitor */
203: eps->nconv = BLZistorr_(blz->istor,"NTEIG",5);
204: EPSMonitor(eps,eps->its,eps->nconv,blz->rstor+BLZistorr_(blz->istor,"IRITZ",5),eps->eigi,blz->rstor+BLZistorr_(blz->istor,"IRITZ",5)+BLZistorr_(blz->istor,"JT",2),BLZistorr_(blz->istor,"NRITZ",5));
205: eps->its = eps->its + 1;
206: if (eps->its >= eps->max_it || eps->nconv >= eps->nev) lflag = 5;
207: break;
208: case 2:
209: /* compute v = B u */
210: for (i=0;i<nvopu;i++) {
211: VecPlaceArray(x,blz->u+i*eps->nloc);
212: VecPlaceArray(y,blz->v+i*eps->nloc);
213: BVApplyMatrix(eps->V,x,y);
214: VecResetArray(x);
215: VecResetArray(y);
216: }
217: break;
218: case 3:
219: /* update shift */
220: PetscInfo1(eps,"Factorization update (sigma=%g)\n",(double)sigma);
221: STSetShift(eps->st,sigma);
222: STGetKSP(eps->st,&ksp);
223: KSPGetPC(ksp,&pc);
224: PCFactorGetMatrix(pc,&M);
225: MatGetInertia(M,&nn,NULL,NULL);
226: PetscBLASIntCast(nn,&nneig);
227: break;
228: case 4:
229: /* copy the initial vector */
230: VecPlaceArray(x,blz->v);
231: BVCopyVec(eps->V,0,x);
232: VecResetArray(x);
233: break;
234: }
236: } while (lflag > 0);
238: BVRestoreArray(eps->V,&pV);
240: eps->nconv = BLZistorr_(blz->istor,"NTEIG",5);
241: eps->reason = EPS_CONVERGED_TOL;
242: if (blz->slice) eps->nev = eps->nconv;
243: for (i=0;i<eps->nconv;i++) eps->eigr[i]=blz->eig[i];
245: if (lflag!=0) {
246: char msg[2048] = "";
247: for (i = 0; i < 33; i++) {
248: if (blz->istor[15] & (1 << i)) { PetscStrcat(msg,blzpack_error[i]); }
249: }
250: SETERRQ2(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15],msg);
251: }
252: VecDestroy(&x);
253: VecDestroy(&y);
254: return(0);
255: }
257: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
258: {
260: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
263: if (!blz->slice && !eps->isgeneralized) {
264: EPSBackTransform_Default(eps);
265: }
266: return(0);
267: }
269: PetscErrorCode EPSReset_BLZPACK(EPS eps)
270: {
272: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
275: PetscFree(blz->istor);
276: PetscFree(blz->rstor);
277: PetscFree(blz->u);
278: PetscFree(blz->v);
279: PetscFree(blz->eig);
280: return(0);
281: }
283: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
284: {
288: PetscFree(eps->data);
289: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",NULL);
290: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackGetBlockSize_C",NULL);
291: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",NULL);
292: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackGetNSteps_C",NULL);
293: return(0);
294: }
296: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
297: {
299: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
300: PetscBool isascii;
303: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
304: if (isascii) {
305: PetscViewerASCIIPrintf(viewer," block size=%d\n",blz->block_size);
306: PetscViewerASCIIPrintf(viewer," maximum number of steps per run=%d\n",blz->nsteps);
307: if (blz->slice) {
308: PetscViewerASCIIPrintf(viewer," computational interval [%f,%f]\n",eps->inta,eps->intb);
309: }
310: }
311: return(0);
312: }
314: PetscErrorCode EPSSetFromOptions_BLZPACK(PetscOptionItems *PetscOptionsObject,EPS eps)
315: {
317: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
318: PetscInt bs,n;
319: PetscBool flg;
322: PetscOptionsHead(PetscOptionsObject,"EPS BLZPACK Options");
324: bs = blz->block_size;
325: PetscOptionsInt("-eps_blzpack_blocksize","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
326: if (flg) { EPSBlzpackSetBlockSize(eps,bs); }
328: n = blz->nsteps;
329: PetscOptionsInt("-eps_blzpack_nsteps","Number of steps","EPSBlzpackSetNSteps",n,&n,&flg);
330: if (flg) { EPSBlzpackSetNSteps(eps,n); }
332: PetscOptionsTail();
333: return(0);
334: }
336: static PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
337: {
339: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
342: if (bs == PETSC_DEFAULT) blz->block_size = 3;
343: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Block size must be positive");
344: else {
345: PetscBLASIntCast(bs,&blz->block_size);
346: }
347: return(0);
348: }
350: /*@
351: EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.
353: Collective on eps
355: Input Parameters:
356: + eps - the eigenproblem solver context
357: - bs - block size
359: Options Database Key:
360: . -eps_blzpack_blocksize - Sets the value of the block size
362: Level: advanced
363: @*/
364: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
365: {
371: PetscTryMethod(eps,"EPSBlzpackSetBlockSize_C",(EPS,PetscInt),(eps,bs));
372: return(0);
373: }
375: static PetscErrorCode EPSBlzpackGetBlockSize_BLZPACK(EPS eps,PetscInt *bs)
376: {
377: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
380: *bs = blz->block_size;
381: return(0);
382: }
384: /*@
385: EPSBlzpackGetBlockSize - Gets the block size used in the BLZPACK solver.
387: Not Collective
389: Input Parameter:
390: . eps - the eigenproblem solver context
392: Output Parameter:
393: . bs - the block size
395: Level: advanced
397: .seealso: EPSBlzpackSetBlockSize()
398: @*/
399: PetscErrorCode EPSBlzpackGetBlockSize(EPS eps,PetscInt *bs)
400: {
406: PetscUseMethod(eps,"EPSBlzpackGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
407: return(0);
408: }
410: static PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
411: {
413: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
416: if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
417: else {
418: PetscBLASIntCast(nsteps,&blz->nsteps);
419: }
420: return(0);
421: }
423: /*@
424: EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
425: package.
427: Collective on eps
429: Input Parameters:
430: + eps - the eigenproblem solver context
431: - nsteps - maximum number of steps
433: Options Database Key:
434: . -eps_blzpack_nsteps - Sets the maximum number of steps per run
436: Level: advanced
438: @*/
439: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
440: {
446: PetscTryMethod(eps,"EPSBlzpackSetNSteps_C",(EPS,PetscInt),(eps,nsteps));
447: return(0);
448: }
450: static PetscErrorCode EPSBlzpackGetNSteps_BLZPACK(EPS eps,PetscInt *nsteps)
451: {
452: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
455: *nsteps = blz->nsteps;
456: return(0);
457: }
459: /*@
460: EPSBlzpackGetNSteps - Gets the maximum number of steps per run in the BLZPACK solver.
462: Not Collective
464: Input Parameter:
465: . eps - the eigenproblem solver context
467: Output Parameter:
468: . nsteps - the maximum number of steps
470: Level: advanced
472: .seealso: EPSBlzpackSetNSteps()
473: @*/
474: PetscErrorCode EPSBlzpackGetNSteps(EPS eps,PetscInt *nsteps)
475: {
481: PetscUseMethod(eps,"EPSBlzpackGetNSteps_C",(EPS,PetscInt*),(eps,nsteps));
482: return(0);
483: }
485: SLEPC_EXTERN PetscErrorCode EPSCreate_BLZPACK(EPS eps)
486: {
488: EPS_BLZPACK *blzpack;
491: PetscNewLog(eps,&blzpack);
492: eps->data = (void*)blzpack;
494: eps->ops->solve = EPSSolve_BLZPACK;
495: eps->ops->setup = EPSSetUp_BLZPACK;
496: eps->ops->setfromoptions = EPSSetFromOptions_BLZPACK;
497: eps->ops->destroy = EPSDestroy_BLZPACK;
498: eps->ops->reset = EPSReset_BLZPACK;
499: eps->ops->view = EPSView_BLZPACK;
500: eps->ops->backtransform = EPSBackTransform_BLZPACK;
502: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",EPSBlzpackSetBlockSize_BLZPACK);
503: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackGetBlockSize_C",EPSBlzpackGetBlockSize_BLZPACK);
504: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",EPSBlzpackSetNSteps_BLZPACK);
505: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackGetNSteps_C",EPSBlzpackGetNSteps_BLZPACK);
506: return(0);
507: }