Actual source code: ex24.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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: */

 11: static char help[] = "Spectrum folding for a standard symmetric eigenproblem.\n\n"
 12:   "The problem matrix is the 2-D Laplacian.\n\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n";

 17: #include <slepceps.h>

 19: /*
 20:    User context for spectrum folding
 21: */
 22: typedef struct {
 23:   Mat       A;
 24:   Vec       w;
 25:   PetscReal target;
 26: } CTX_FOLD;

 28: /*
 29:    Auxiliary routines
 30: */
 31: PetscErrorCode MatMult_Fold(Mat,Vec,Vec);
 32: PetscErrorCode RayleighQuotient(Mat,Vec,PetscScalar*);
 33: PetscErrorCode ComputeResidualNorm(Mat,PetscScalar,Vec,PetscReal*);

 35: int main(int argc,char **argv)
 36: {
 37:   Mat            A,M,P;       /* problem matrix, shell matrix and preconditioner */
 38:   Vec            x;           /* eigenvector */
 39:   EPS            eps;         /* eigenproblem solver context */
 40:   ST             st;          /* spectral transformation */
 41:   KSP            ksp;
 42:   PC             pc;
 43:   EPSType        type;
 44:   CTX_FOLD       *ctx;
 45:   PetscInt       nconv,N,n=10,m,nloc,mloc,Istart,Iend,II,i,j;
 46:   PetscReal      *error,*evals,target=0.0,tol;
 47:   PetscScalar    lambda;
 48:   PetscBool      flag,terse,errok,hasmat;

 50:   SlepcInitialize(&argc,&argv,(char*)0,help);

 52:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 53:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 54:   if (!flag) m=n;
 55:   PetscOptionsGetReal(NULL,NULL,"-target",&target,NULL);
 56:   N = n*m;
 57:   PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum Folding, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid) target=%f\n\n",N,n,m,(double)target);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:      Compute the 5-point stencil Laplacian
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   MatCreate(PETSC_COMM_WORLD,&A);
 64:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 65:   MatSetFromOptions(A);
 66:   MatSetUp(A);

 68:   MatGetOwnershipRange(A,&Istart,&Iend);
 69:   for (II=Istart;II<Iend;II++) {
 70:     i = II/n; j = II-i*n;
 71:     if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
 72:     if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
 73:     if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
 74:     if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
 75:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 76:   }

 78:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 80:   MatCreateVecs(A,&x,NULL);
 81:   MatGetLocalSize(A,&nloc,&mloc);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:                 Create shell matrix to perform spectrum folding
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 86:   PetscNew(&ctx);
 87:   ctx->A = A;
 88:   ctx->target = target;
 89:   VecDuplicate(x,&ctx->w);

 91:   MatCreateShell(PETSC_COMM_WORLD,nloc,mloc,N,N,ctx,&M);
 92:   MatShellSetOperation(M,MATOP_MULT,(void(*)(void))MatMult_Fold);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:                 Create the eigensolver and set various options
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   EPSCreate(PETSC_COMM_WORLD,&eps);
 99:   EPSSetOperators(eps,M,NULL);
100:   EPSSetProblemType(eps,EPS_HEP);
101:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
102:   EPSSetFromOptions(eps);

104:   PetscObjectTypeCompareAny((PetscObject)eps,&flag,EPSGD,EPSJD,EPSBLOPEX,EPSLOBPCG,EPSRQCG,"");
105:   if (flag) {
106:     /*
107:        Build preconditioner specific for this application (diagonal of A^2)
108:     */
109:     MatGetRowSum(A,x);
110:     VecScale(x,-1.0);
111:     VecShift(x,20.0);
112:     MatCreate(PETSC_COMM_WORLD,&P);
113:     MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,N,N);
114:     MatSetFromOptions(P);
115:     MatSetUp(P);
116:     MatDiagonalSet(P,x,INSERT_VALUES);
117:     /*
118:        Set diagonal preconditioner
119:     */
120:     EPSGetST(eps,&st);
121:     STSetType(st,STPRECOND);
122:     STSetPreconditionerMat(st,P);
123:     MatDestroy(&P);
124:     STGetKSP(st,&ksp);
125:     KSPGetPC(ksp,&pc);
126:     PCSetType(pc,PCJACOBI);
127:     STPrecondGetKSPHasMat(st,&hasmat);
128:     PetscPrintf(PETSC_COMM_WORLD," Preconditioned solver, hasmat=%s\n",hasmat?"true":"false");
129:   }

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:                       Solve the eigensystem
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

135:   EPSSolve(eps);
136:   EPSGetType(eps,&type);
137:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
138:   EPSGetTolerances(eps,&tol,NULL);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:                     Display solution and clean up
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

144:   EPSGetConverged(eps,&nconv);
145:   PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv);
146:   if (nconv>0) {
147:     PetscMalloc2(nconv,&evals,nconv,&error);
148:     for (i=0;i<nconv;i++) {
149:       /*  Get i-th eigenvector, compute eigenvalue approximation from
150:           Rayleigh quotient and compute residual norm */
151:       EPSGetEigenpair(eps,i,NULL,NULL,x,NULL);
152:       RayleighQuotient(A,x,&lambda);
153:       ComputeResidualNorm(A,lambda,x,&error[i]);
154: #if defined(PETSC_USE_COMPLEX)
155:       evals[i] = PetscRealPart(lambda);
156: #else
157:       evals[i] = lambda;
158: #endif
159:     }
160:     PetscOptionsHasName(NULL,NULL,"-terse",&terse);
161:     if (!terse) {
162:       PetscCall(PetscPrintf(PETSC_COMM_WORLD,
163:            "           k              ||Ax-kx||\n"
164:            "   ----------------- ------------------\n"));
165:       for (i=0;i<nconv;i++) PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12.2g\n",(double)evals[i],(double)error[i]);
166:     } else {
167:       errok = PETSC_TRUE;
168:       for (i=0;i<nconv;i++) errok = (errok && error[i]<5.0*tol)? PETSC_TRUE: PETSC_FALSE;
169:       if (!errok) PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nconv);
170:       else {
171:         PetscPrintf(PETSC_COMM_WORLD," nconv=%" PetscInt_FMT " eigenvalues computed up to the required tolerance:",nconv);
172:         for (i=0;i<nconv;i++) PetscPrintf(PETSC_COMM_WORLD," %.5f",(double)evals[i]);
173:       }
174:     }
175:     PetscPrintf(PETSC_COMM_WORLD,"\n");
176:     PetscFree2(evals,error);
177:   }

179:   EPSDestroy(&eps);
180:   MatDestroy(&A);
181:   MatDestroy(&M);
182:   VecDestroy(&ctx->w);
183:   VecDestroy(&x);
184:   PetscFree(ctx);
185:   SlepcFinalize();
186:   return 0;
187: }

189: /*
190:     Matrix-vector product subroutine for the spectrum folding.
191:        y <-- (A-t*I)^2*x
192:  */
193: PetscErrorCode MatMult_Fold(Mat M,Vec x,Vec y)
194: {
195:   CTX_FOLD       *ctx;
196:   PetscScalar    sigma;

199:   MatShellGetContext(M,&ctx);
200:   sigma = -ctx->target;
201:   MatMult(ctx->A,x,ctx->w);
202:   VecAXPY(ctx->w,sigma,x);
203:   MatMult(ctx->A,ctx->w,y);
204:   VecAXPY(y,sigma,ctx->w);
205:   PetscFunctionReturn(0);
206: }

208: /*
209:     Computes the Rayleigh quotient of a vector x
210:        r <-- x^T*A*x       (assumes x has unit norm)
211:  */
212: PetscErrorCode RayleighQuotient(Mat A,Vec x,PetscScalar *r)
213: {
214:   Vec            Ax;

217:   VecDuplicate(x,&Ax);
218:   MatMult(A,x,Ax);
219:   VecDot(Ax,x,r);
220:   VecDestroy(&Ax);
221:   PetscFunctionReturn(0);
222: }

224: /*
225:     Computes the residual norm of an approximate eigenvector x, |A*x-lambda*x|
226:  */
227: PetscErrorCode ComputeResidualNorm(Mat A,PetscScalar lambda,Vec x,PetscReal *r)
228: {
229:   Vec            Ax;

232:   VecDuplicate(x,&Ax);
233:   MatMult(A,x,Ax);
234:   VecAXPY(Ax,-lambda,x);
235:   VecNorm(Ax,NORM_2,r);
236:   VecDestroy(&Ax);
237:   PetscFunctionReturn(0);
238: }

240: /*TEST

242:    testset:
243:       args: -n 15 -eps_nev 1 -eps_ncv 12 -eps_max_it 1000 -eps_tol 1e-5 -terse
244:       filter: grep -v Solution
245:       test:
246:          suffix: 1
247:       test:
248:          suffix: 1_lobpcg
249:          args: -eps_type lobpcg
250:          requires: !single
251:       test:
252:          suffix: 1_gd
253:          args: -eps_type gd
254:          requires: !single

256: TEST*/