Actual source code: test39.c
slepc-3.19.0 2023-03-31
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;
26: PetscFunctionBeginUser;
27: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
28: nloc = PETSC_DECIDE;
29: PetscCall(PetscSplitOwnership(comm,&nloc,&N));
30: if (rank==0) nloc += shift;
31: else if (rank==1) nloc -= shift;
33: PetscCall(MatCreate(comm,A));
34: PetscCall(MatSetSizes(*A,nloc,nloc,N,N));
35: PetscCall(MatSetFromOptions(*A));
36: PetscCall(MatSetUp(*A));
37: PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
38: for (II=Istart;II<Iend;II++) {
39: i = II/n; j = II-i*n;
40: if (i>0) PetscCall(MatSetValue(*A,II,II-n,-1.0,INSERT_VALUES));
41: if (i<m-1) PetscCall(MatSetValue(*A,II,II+n,-1.0,INSERT_VALUES));
42: if (j>0) PetscCall(MatSetValue(*A,II,II-1,-1.0,INSERT_VALUES));
43: if (j<n-1) PetscCall(MatSetValue(*A,II,II+1,-1.0,INSERT_VALUES));
44: PetscCall(MatSetValue(*A,II,II,4.0,INSERT_VALUES));
45: }
46: PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
47: PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
48: PetscFunctionReturn(PETSC_SUCCESS);
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: PetscFunctionBeginUser;
60: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
61: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
62: PetscCheck(size>1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least two processes");
63: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
64: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
65: N = n*m;
66: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Create 2-D Laplacian matrices
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,1,&A));
73: PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,-1,&B));
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Create the eigensolver, set options and solve the eigensystem
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First solve:\n\n"));
80: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
81: PetscCall(EPSSetOperators(eps,A,NULL));
82: PetscCall(EPSSetProblemType(eps,EPS_HEP));
83: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
84: PetscCall(EPSSetDimensions(eps,nev,PETSC_DEFAULT,PETSC_DEFAULT));
85: PetscCall(EPSSetFromOptions(eps));
87: PetscCall(EPSSolve(eps));
89: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90: Display solution of first solve
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
94: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
95: else {
96: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
97: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
98: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
99: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
100: }
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Solve with second matrix
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSecond solve:\n\n"));
107: /* PetscCall(EPSReset(eps)); */ /* not required, will be called in EPSSetOperators() */
108: PetscCall(EPSSetOperators(eps,B,NULL));
109: PetscCall(EPSSolve(eps));
111: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
112: else {
113: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
114: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
115: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
116: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
117: }
119: PetscCall(EPSDestroy(&eps));
120: PetscCall(MatDestroy(&A));
121: PetscCall(MatDestroy(&B));
122: PetscCall(SlepcFinalize());
123: return 0;
124: }
126: /*TEST
128: testset:
129: nsize: 2
130: requires: !single
131: output_file: output/test39_1.out
132: test:
133: suffix: 1
134: args: -eps_type {{krylovschur arnoldi lobpcg lapack}} -terse
135: test:
136: suffix: 1_lanczos
137: args: -eps_type lanczos -eps_lanczos_reorthog local -terse
139: TEST*/