Actual source code: test4.c
slepc-3.18.1 2022-11-02
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;
31: SlepcInitialize(&argc,&argv,(char*)0,help);
32: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
33: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
34: if (!flag) m=n;
35: N = n*m;
36: 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: MatCreate(PETSC_COMM_WORLD,&K);
44: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
45: MatSetFromOptions(K);
46: MatSetUp(K);
47: MatGetOwnershipRange(K,&Istart,&Iend);
48: for (II=Istart;II<Iend;II++) {
49: i = II/n; j = II-i*n;
50: if (i>0) MatSetValue(K,II,II-n,-1.0,INSERT_VALUES);
51: if (i<m-1) MatSetValue(K,II,II+n,-1.0,INSERT_VALUES);
52: if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
53: if (j<n-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
54: MatSetValue(K,II,II,4.0,INSERT_VALUES);
55: }
56: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
59: /* C is the 1-D Laplacian on horizontal lines */
60: MatCreate(PETSC_COMM_WORLD,&C);
61: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
62: MatSetFromOptions(C);
63: MatSetUp(C);
64: MatGetOwnershipRange(C,&Istart,&Iend);
65: for (II=Istart;II<Iend;II++) {
66: i = II/n; j = II-i*n;
67: if (j>0) MatSetValue(C,II,II-1,-1.0,INSERT_VALUES);
68: if (j<n-1) MatSetValue(C,II,II+1,-1.0,INSERT_VALUES);
69: MatSetValue(C,II,II,2.0,INSERT_VALUES);
70: }
71: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
74: /* M is a diagonal matrix */
75: MatCreate(PETSC_COMM_WORLD,&M);
76: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
77: MatSetFromOptions(M);
78: MatSetUp(M);
79: MatGetOwnershipRange(M,&Istart,&Iend);
80: for (II=Istart;II<Iend;II++) MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
81: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create a standalone EPS with appropriate settings
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: EPSCreate(PETSC_COMM_WORLD,&eps);
89: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
90: #if defined(PETSC_USE_COMPLEX)
91: EPSSetTarget(eps,0.01*PETSC_i);
92: #endif
93: EPSGetST(eps,&st);
94: STSetType(st,STSINVERT);
95: STGetKSP(st,&ksp);
96: KSPSetType(ksp,KSPBCGS);
97: KSPGetPC(ksp,&pc);
98: PCSetType(pc,PCJACOBI);
99: EPSSetFromOptions(eps);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Create the eigensolver and solve the eigensystem
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: PEPCreate(PETSC_COMM_WORLD,&pep);
106: PetscObjectSetName((PetscObject)pep,"PEP_solver");
107: A[0] = K; A[1] = C; A[2] = M;
108: PEPSetOperators(pep,3,A);
109: PEPSetType(pep,PEPLINEAR);
110: PEPSetProblemType(pep,PEP_GENERAL);
111: PEPLinearSetEPS(pep,eps);
112: PEPSetFromOptions(pep);
113: PEPSolve(pep);
114: PEPLinearGetLinearization(pep,&alpha,&beta);
115: PetscPrintf(PETSC_COMM_WORLD," Linearization with alpha=%g, beta=%g",(double)alpha,(double)beta);
116: PEPLinearGetExplicitMatrix(pep,&expmat);
117: if (expmat) PetscPrintf(PETSC_COMM_WORLD," with explicit matrix");
118: PetscPrintf(PETSC_COMM_WORLD,"\n");
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Display solution and clean up
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
125: PEPDestroy(&pep);
126: EPSDestroy(&eps);
127: MatDestroy(&M);
128: MatDestroy(&C);
129: MatDestroy(&K);
130: 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*/