Actual source code: test4.c
slepc-3.19.2 2023-09-05
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[] = "Solve a quadratic problem with PEPLINEAR with a user-provided EPS.\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 <slepcpep.h>
18: int main(int argc,char **argv)
19: {
20: Mat M,C,K,A[3];
21: PEP pep;
22: PetscInt N,n=10,m,Istart,Iend,II,i,j;
23: PetscBool flag,expmat;
24: PetscReal alpha,beta;
25: EPS eps;
26: ST st;
27: KSP ksp;
28: PC pc;
30: PetscFunctionBeginUser;
31: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
32: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
33: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
34: if (!flag) m=n;
35: N = n*m;
36: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
38: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
42: /* K is the 2-D Laplacian */
43: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
44: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
45: PetscCall(MatSetFromOptions(K));
46: PetscCall(MatSetUp(K));
47: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
48: for (II=Istart;II<Iend;II++) {
49: i = II/n; j = II-i*n;
50: if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
51: if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
52: if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
53: if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
54: PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
55: }
56: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
59: /* C is the 1-D Laplacian on horizontal lines */
60: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
61: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
62: PetscCall(MatSetFromOptions(C));
63: PetscCall(MatSetUp(C));
64: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
65: for (II=Istart;II<Iend;II++) {
66: i = II/n; j = II-i*n;
67: if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
68: if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
69: PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
70: }
71: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
72: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
74: /* M is a diagonal matrix */
75: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
76: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
77: PetscCall(MatSetFromOptions(M));
78: PetscCall(MatSetUp(M));
79: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
80: for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
81: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
82: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create a standalone EPS with appropriate settings
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
89: PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
90: #if defined(PETSC_USE_COMPLEX)
91: PetscCall(EPSSetTarget(eps,0.01*PETSC_i));
92: #endif
93: PetscCall(EPSGetST(eps,&st));
94: PetscCall(STSetType(st,STSINVERT));
95: PetscCall(STGetKSP(st,&ksp));
96: PetscCall(KSPSetType(ksp,KSPBCGS));
97: PetscCall(KSPGetPC(ksp,&pc));
98: PetscCall(PCSetType(pc,PCJACOBI));
99: PetscCall(EPSSetFromOptions(eps));
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Create the eigensolver and solve the eigensystem
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
106: PetscCall(PetscObjectSetName((PetscObject)pep,"PEP_solver"));
107: A[0] = K; A[1] = C; A[2] = M;
108: PetscCall(PEPSetOperators(pep,3,A));
109: PetscCall(PEPSetType(pep,PEPLINEAR));
110: PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
111: PetscCall(PEPLinearSetEPS(pep,eps));
112: PetscCall(PEPSetFromOptions(pep));
113: PetscCall(PEPSolve(pep));
114: PetscCall(PEPLinearGetLinearization(pep,&alpha,&beta));
115: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Linearization with alpha=%g, beta=%g",(double)alpha,(double)beta));
116: PetscCall(PEPLinearGetExplicitMatrix(pep,&expmat));
117: if (expmat) PetscCall(PetscPrintf(PETSC_COMM_WORLD," with explicit matrix"));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Display solution and clean up
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
125: PetscCall(PEPDestroy(&pep));
126: PetscCall(EPSDestroy(&eps));
127: PetscCall(MatDestroy(&M));
128: PetscCall(MatDestroy(&C));
129: PetscCall(MatDestroy(&K));
130: PetscCall(SlepcFinalize());
131: return 0;
132: }
134: /*TEST
136: testset:
137: args: -pep_linear_explicitmatrix -pep_view_vectors ::ascii_info
138: test:
139: suffix: 1_real
140: requires: !single !complex
141: test:
142: suffix: 1
143: requires: !single complex
145: TEST*/