Actual source code: test39.c
slepc-3.22.1 2024-10-28
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(MatGetOwnershipRange(*A,&Istart,&Iend));
37: for (II=Istart;II<Iend;II++) {
38: i = II/n; j = II-i*n;
39: if (i>0) PetscCall(MatSetValue(*A,II,II-n,-1.0,INSERT_VALUES));
40: if (i<m-1) PetscCall(MatSetValue(*A,II,II+n,-1.0,INSERT_VALUES));
41: if (j>0) PetscCall(MatSetValue(*A,II,II-1,-1.0,INSERT_VALUES));
42: if (j<n-1) PetscCall(MatSetValue(*A,II,II+1,-1.0,INSERT_VALUES));
43: PetscCall(MatSetValue(*A,II,II,4.0,INSERT_VALUES));
44: }
45: PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
46: PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
47: PetscFunctionReturn(PETSC_SUCCESS);
48: }
50: int main(int argc,char **argv)
51: {
52: Mat A,B;
53: EPS eps;
54: PetscInt N,n=10,m=11,nev=3;
55: PetscMPIInt size;
56: PetscBool flag,terse;
58: PetscFunctionBeginUser;
59: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
60: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
61: PetscCheck(size>1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least two processes");
62: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
63: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
64: N = n*m;
65: PetscCall(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: PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,1,&A));
72: PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,-1,&B));
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create the eigensolver, set options and solve the eigensystem
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First solve:\n\n"));
79: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
80: PetscCall(EPSSetOperators(eps,A,NULL));
81: PetscCall(EPSSetProblemType(eps,EPS_HEP));
82: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
83: PetscCall(EPSSetDimensions(eps,nev,PETSC_DETERMINE,PETSC_DETERMINE));
84: PetscCall(EPSSetFromOptions(eps));
86: PetscCall(EPSSolve(eps));
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Display solution of first solve
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
93: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
94: else {
95: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
96: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
97: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
98: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
99: }
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Solve with second matrix
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSecond solve:\n\n"));
106: /* PetscCall(EPSReset(eps)); */ /* not required, will be called in EPSSetOperators() */
107: PetscCall(EPSSetOperators(eps,B,NULL));
108: PetscCall(EPSSolve(eps));
110: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
111: else {
112: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
113: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
114: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
115: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
116: }
118: PetscCall(EPSDestroy(&eps));
119: PetscCall(MatDestroy(&A));
120: PetscCall(MatDestroy(&B));
121: PetscCall(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*/