Actual source code: test39.c

slepc-3.17.0 2022-03-31
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[] = "Tests multiple calls to EPSSolve with matrices of different local size.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepceps.h>

 18: /*
 19:    Create 2-D Laplacian matrix
 20: */
 21: PetscErrorCode Laplacian(MPI_Comm comm,PetscInt n,PetscInt m,PetscInt shift,Mat *A)
 22: {
 23:   PetscInt       N = n*m,i,j,II,Istart,Iend,nloc;
 24:   PetscMPIInt    rank;

 27:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 28:   nloc = PETSC_DECIDE;
 29:   PetscSplitOwnership(comm,&nloc,&N);
 30:   if (rank==0) nloc += shift;
 31:   else if (rank==1) nloc -= shift;

 33:   MatCreate(comm,A);
 34:   MatSetSizes(*A,nloc,nloc,N,N);
 35:   MatSetFromOptions(*A);
 36:   MatSetUp(*A);
 37:   MatGetOwnershipRange(*A,&Istart,&Iend);
 38:   for (II=Istart;II<Iend;II++) {
 39:     i = II/n; j = II-i*n;
 40:     if (i>0) MatSetValue(*A,II,II-n,-1.0,INSERT_VALUES);
 41:     if (i<m-1) MatSetValue(*A,II,II+n,-1.0,INSERT_VALUES);
 42:     if (j>0) MatSetValue(*A,II,II-1,-1.0,INSERT_VALUES);
 43:     if (j<n-1) MatSetValue(*A,II,II+1,-1.0,INSERT_VALUES);
 44:     MatSetValue(*A,II,II,4.0,INSERT_VALUES);
 45:   }
 46:   MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
 47:   MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
 48:   PetscFunctionReturn(0);
 49: }

 51: int main(int argc,char **argv)
 52: {
 53:   Mat            A,B;
 54:   EPS            eps;
 55:   PetscInt       N,n=10,m=11,nev=3;
 56:   PetscMPIInt    size;
 57:   PetscBool      flag,terse;

 59:   SlepcInitialize(&argc,&argv,(char*)0,help);
 60:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 62:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 63:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 64:   N = n*m;
 65:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:                 Create 2-D Laplacian matrices
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 71:   Laplacian(PETSC_COMM_WORLD,n,m,1,&A);
 72:   Laplacian(PETSC_COMM_WORLD,n,m,-1,&B);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:         Create the eigensolver, set options and solve the eigensystem
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   PetscPrintf(PETSC_COMM_WORLD,"First solve:\n\n");
 79:   EPSCreate(PETSC_COMM_WORLD,&eps);
 80:   EPSSetOperators(eps,A,NULL);
 81:   EPSSetProblemType(eps,EPS_HEP);
 82:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
 83:   EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT);
 84:   EPSSetFromOptions(eps);

 86:   EPSSolve(eps);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                     Display solution of first solve
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
 93:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 94:   else {
 95:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 96:     EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
 97:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
 98:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 99:   }

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                        Solve with second matrix
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   PetscPrintf(PETSC_COMM_WORLD,"\nSecond solve:\n\n");
106:   /*EPSReset(eps);*/  /* not required, will be called in EPSSetOperators() */
107:   EPSSetOperators(eps,B,NULL);
108:   EPSSolve(eps);

110:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
111:   else {
112:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
113:     EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
114:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
115:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
116:   }

118:   EPSDestroy(&eps);
119:   MatDestroy(&A);
120:   MatDestroy(&B);
121:   SlepcFinalize();
122:   return 0;
123: }

125: /*TEST

127:    testset:
128:       nsize: 2
129:       requires: !single
130:       output_file: output/test39_1.out
131:       test:
132:          suffix: 1
133:          args: -eps_type {{krylovschur arnoldi lobpcg lapack}} -terse
134:       test:
135:          suffix: 1_lanczos
136:          args: -eps_type lanczos -eps_lanczos_reorthog local -terse

138: TEST*/