Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
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";
15 :
16 : #include <slepceps.h>
17 :
18 : /*
19 : Create 2-D Laplacian matrix
20 : */
21 20 : PetscErrorCode Laplacian(MPI_Comm comm,PetscInt n,PetscInt m,PetscInt shift,Mat *A)
22 : {
23 20 : PetscInt N = n*m,i,j,II,Istart,Iend,nloc;
24 20 : PetscMPIInt rank;
25 :
26 20 : PetscFunctionBeginUser;
27 20 : PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
28 20 : nloc = PETSC_DECIDE;
29 20 : PetscCall(PetscSplitOwnership(comm,&nloc,&N));
30 20 : if (rank==0) nloc += shift;
31 10 : else if (rank==1) nloc -= shift;
32 :
33 20 : PetscCall(MatCreate(comm,A));
34 20 : PetscCall(MatSetSizes(*A,nloc,nloc,N,N));
35 20 : PetscCall(MatSetFromOptions(*A));
36 20 : PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
37 1120 : for (II=Istart;II<Iend;II++) {
38 1100 : i = II/n; j = II-i*n;
39 1100 : if (i>0) PetscCall(MatSetValue(*A,II,II-n,-1.0,INSERT_VALUES));
40 1100 : if (i<m-1) PetscCall(MatSetValue(*A,II,II+n,-1.0,INSERT_VALUES));
41 1100 : if (j>0) PetscCall(MatSetValue(*A,II,II-1,-1.0,INSERT_VALUES));
42 1100 : if (j<n-1) PetscCall(MatSetValue(*A,II,II+1,-1.0,INSERT_VALUES));
43 1100 : PetscCall(MatSetValue(*A,II,II,4.0,INSERT_VALUES));
44 : }
45 20 : PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
46 20 : PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
47 20 : PetscFunctionReturn(PETSC_SUCCESS);
48 : }
49 :
50 10 : int main(int argc,char **argv)
51 : {
52 10 : Mat A,B;
53 10 : EPS eps;
54 10 : PetscInt N,n=10,m=11,nev=3;
55 10 : PetscMPIInt size;
56 10 : PetscBool flag,terse;
57 :
58 10 : PetscFunctionBeginUser;
59 10 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
60 10 : PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
61 10 : PetscCheck(size>1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least two processes");
62 10 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
63 10 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
64 10 : N = n*m;
65 10 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
66 :
67 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68 : Create 2-D Laplacian matrices
69 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70 :
71 10 : PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,1,&A));
72 10 : PetscCall(Laplacian(PETSC_COMM_WORLD,n,m,-1,&B));
73 :
74 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75 : Create the eigensolver, set options and solve the eigensystem
76 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77 :
78 10 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"First solve:\n\n"));
79 10 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
80 10 : PetscCall(EPSSetOperators(eps,A,NULL));
81 10 : PetscCall(EPSSetProblemType(eps,EPS_HEP));
82 10 : PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
83 10 : PetscCall(EPSSetDimensions(eps,nev,PETSC_DETERMINE,PETSC_DETERMINE));
84 10 : PetscCall(EPSSetFromOptions(eps));
85 :
86 10 : PetscCall(EPSSolve(eps));
87 :
88 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89 : Display solution of first solve
90 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91 :
92 10 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
93 10 : if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
94 : else {
95 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
96 0 : PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
97 0 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
98 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
99 : }
100 :
101 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102 : Solve with second matrix
103 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104 :
105 10 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSecond solve:\n\n"));
106 : /* PetscCall(EPSReset(eps)); */ /* not required, will be called in EPSSetOperators() */
107 10 : PetscCall(EPSSetOperators(eps,B,NULL));
108 10 : PetscCall(EPSSolve(eps));
109 :
110 10 : if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
111 : else {
112 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
113 0 : PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
114 0 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
115 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
116 : }
117 :
118 10 : PetscCall(EPSDestroy(&eps));
119 10 : PetscCall(MatDestroy(&A));
120 10 : PetscCall(MatDestroy(&B));
121 10 : PetscCall(SlepcFinalize());
122 : return 0;
123 : }
124 :
125 : /*TEST
126 :
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
137 :
138 : TEST*/
|