Actual source code: ex19.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[] = "Standard symmetric eigenproblem for the 3-D Laplacian built with the DM interface.\n\n"
 12: "Use -seed <k> to modify the random initial vector.\n"
 13: "Use -da_grid_x <nx> etc. to change the problem size.\n\n";

 15: #include <slepceps.h>
 16: #include <petscdmda.h>
 17: #include <petsctime.h>

 19: PetscErrorCode GetExactEigenvalues(PetscInt M,PetscInt N,PetscInt P,PetscInt nconv,PetscReal *exact)
 20: {
 21:   PetscInt       n,i,j,k,l;
 22:   PetscReal      *evals,ax,ay,az,sx,sy,sz;

 25:   ax = PETSC_PI/2/(M+1);
 26:   ay = PETSC_PI/2/(N+1);
 27:   az = PETSC_PI/2/(P+1);
 28:   n = PetscCeilReal(PetscPowReal((PetscReal)nconv,0.33333)+1);
 29:   PetscMalloc1(n*n*n,&evals);
 30:   l = 0;
 31:   for (i=1;i<=n;i++) {
 32:     sx = PetscSinReal(ax*i);
 33:     for (j=1;j<=n;j++) {
 34:       sy = PetscSinReal(ay*j);
 35:       for (k=1;k<=n;k++) {
 36:         sz = PetscSinReal(az*k);
 37:         evals[l++] = 4.0*(sx*sx+sy*sy+sz*sz);
 38:       }
 39:     }
 40:   }
 41:   PetscSortReal(n*n*n,evals);
 42:   for (i=0;i<nconv;i++) exact[i] = evals[i];
 43:   PetscFree(evals);
 44:   PetscFunctionReturn(0);
 45: }

 47: PetscErrorCode FillMatrix(DM da,Mat A)
 48: {
 49:   PetscInt       i,j,k,mx,my,mz,xm,ym,zm,xs,ys,zs,idx;
 50:   PetscScalar    v[7];
 51:   MatStencil     row,col[7];

 54:   DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);
 55:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

 57:   for (k=zs;k<zs+zm;k++) {
 58:     for (j=ys;j<ys+ym;j++) {
 59:       for (i=xs;i<xs+xm;i++) {
 60:         row.i=i; row.j=j; row.k=k;
 61:         col[0].i=row.i; col[0].j=row.j; col[0].k=row.k;
 62:         v[0]=6.0;
 63:         idx=1;
 64:         if (k>0) { v[idx]=-1.0; col[idx].i=i; col[idx].j=j; col[idx].k=k-1; idx++; }
 65:         if (j>0) { v[idx]=-1.0; col[idx].i=i; col[idx].j=j-1; col[idx].k=k; idx++; }
 66:         if (i>0) { v[idx]=-1.0; col[idx].i=i-1; col[idx].j=j; col[idx].k=k; idx++; }
 67:         if (i<mx-1) { v[idx]=-1.0; col[idx].i=i+1; col[idx].j=j; col[idx].k=k; idx++; }
 68:         if (j<my-1) { v[idx]=-1.0; col[idx].i=i; col[idx].j=j+1; col[idx].k=k; idx++; }
 69:         if (k<mz-1) { v[idx]=-1.0; col[idx].i=i; col[idx].j=j; col[idx].k=k+1; idx++; }
 70:         MatSetValuesStencil(A,1,&row,idx,col,v,INSERT_VALUES);
 71:       }
 72:     }
 73:   }
 74:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 75:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 76:   PetscFunctionReturn(0);
 77: }

 79: int main(int argc,char **argv)
 80: {
 81:   Mat            A;               /* operator matrix */
 82:   EPS            eps;             /* eigenproblem solver context */
 83:   EPSType        type;
 84:   DM             da;
 85:   Vec            v0;
 86:   PetscReal      error,tol,re,im,*exact;
 87:   PetscScalar    kr,ki;
 88:   PetscInt       M,N,P,m,n,p,nev,maxit,i,its,nconv,seed;
 89:   PetscLogDouble t1,t2,t3;
 90:   PetscBool      flg,terse;
 91:   PetscRandom    rctx;

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

 95:   PetscPrintf(PETSC_COMM_WORLD,"\n3-D Laplacian Eigenproblem\n\n");

 97:   /* show detailed info unless -terse option is given by user */
 98:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Compute the operator matrix that defines the eigensystem, Ax=kx
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

104:   PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,
105:                       DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,10,10,10,
106:                       PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,
107:                       1,1,NULL,NULL,NULL,&da));
108:   DMSetFromOptions(da);
109:   DMSetUp(da);

111:   /* print DM information */
112:   DMDAGetInfo(da,NULL,&M,&N,&P,&m,&n,&p,NULL,NULL,NULL,NULL,NULL,NULL);
113:   PetscPrintf(PETSC_COMM_WORLD," Grid partitioning: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",m,n,p);

115:   /* create and fill the matrix */
116:   DMCreateMatrix(da,&A);
117:   FillMatrix(da,A);

119:   /* create random initial vector */
120:   seed = 1;
121:   PetscOptionsGetInt(NULL,NULL,"-seed",&seed,NULL);
123:   MatCreateVecs(A,&v0,NULL);
124:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
125:   PetscRandomSetFromOptions(rctx);
126:   for (i=0;i<seed;i++) {   /* simulate different seeds in the random generator */
127:     VecSetRandom(v0,rctx);
128:   }

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:                 Create the eigensolver and set various options
132:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

134:   /*
135:      Create eigensolver context
136:   */
137:   EPSCreate(PETSC_COMM_WORLD,&eps);

139:   /*
140:      Set operators. In this case, it is a standard eigenvalue problem
141:   */
142:   EPSSetOperators(eps,A,NULL);
143:   EPSSetProblemType(eps,EPS_HEP);

145:   /*
146:      Set specific solver options
147:   */
148:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
149:   EPSSetTolerances(eps,1e-8,PETSC_DEFAULT);
150:   EPSSetInitialSpace(eps,1,&v0);

152:   /*
153:      Set solver parameters at runtime
154:   */
155:   EPSSetFromOptions(eps);

157:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158:                       Solve the eigensystem
159:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

161:   PetscTime(&t1);
162:   EPSSetUp(eps);
163:   PetscTime(&t2);
164:   EPSSolve(eps);
165:   PetscTime(&t3);
166:   if (!terse) {
167:     EPSGetIterationNumber(eps,&its);
168:     PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its);

170:     /*
171:        Optional: Get some information from the solver and display it
172:     */
173:     EPSGetType(eps,&type);
174:     PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
175:     EPSGetDimensions(eps,&nev,NULL,NULL);
176:     PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
177:     EPSGetTolerances(eps,&tol,&maxit);
178:     PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit);
179:   }

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:                     Display solution and clean up
183:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

185:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
186:   else {
187:     /*
188:        Get number of converged approximate eigenpairs
189:     */
190:     EPSGetConverged(eps,&nconv);
191:     PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %" PetscInt_FMT "\n\n",nconv);

193:     if (nconv>0) {
194:       PetscMalloc1(nconv,&exact);
195:       GetExactEigenvalues(M,N,P,nconv,exact);
196:       /*
197:          Display eigenvalues and relative errors
198:       */
199:       PetscCall(PetscPrintf(PETSC_COMM_WORLD,
200:            "           k          ||Ax-kx||/||kx||   Eigenvalue Error \n"
201:            "   ----------------- ------------------ ------------------\n"));

203:       for (i=0;i<nconv;i++) {
204:         /*
205:           Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
206:           ki (imaginary part)
207:         */
208:         EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL);
209:         /*
210:            Compute the relative error associated to each eigenpair
211:         */
212:         EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);

214: #if defined(PETSC_USE_COMPLEX)
215:         re = PetscRealPart(kr);
216:         im = PetscImaginaryPart(kr);
217: #else
218:         re = kr;
219:         im = ki;
220: #endif
222:         PetscPrintf(PETSC_COMM_WORLD,"   %12g       %12g        %12g\n",(double)re,(double)error,(double)PetscAbsReal(re-exact[i]));
223:       }
224:       PetscFree(exact);
225:       PetscPrintf(PETSC_COMM_WORLD,"\n");
226:     }
227:   }

229:   /*
230:      Show computing times
231:   */
232:   PetscOptionsHasName(NULL,NULL,"-showtimes",&flg);
233:   if (flg) PetscPrintf(PETSC_COMM_WORLD," Elapsed time: %g (setup), %g (solve)\n",(double)(t2-t1),(double)(t3-t2));

235:   /*
236:      Free work space
237:   */
238:   EPSDestroy(&eps);
239:   MatDestroy(&A);
240:   VecDestroy(&v0);
241:   PetscRandomDestroy(&rctx);
242:   DMDestroy(&da);
243:   SlepcFinalize();
244:   return 0;
245: }

247: /*TEST

249:    testset:
250:       args: -eps_nev 8 -terse
251:       requires: double
252:       output_file: output/ex19_1.out
253:       test:
254:          suffix: 1_krylovschur
255:          args: -eps_type krylovschur -eps_ncv 64
256:       test:
257:          suffix: 1_lobpcg
258:          args: -eps_type lobpcg -eps_tol 1e-7
259:       test:
260:          suffix: 1_blopex
261:          args: -eps_type blopex -eps_tol 1e-7 -eps_blopex_blocksize 4
262:          requires: blopex

264: TEST*/