Actual source code: test4.c
slepc-3.9.1 2018-05-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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,cform;
23: PetscBool flag,expmat;
24: EPS eps;
25: ST st;
26: KSP ksp;
27: PC pc;
30: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
32: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
33: if (!flag) m=n;
34: N = n*m;
35: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
37: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41: /* K is the 2-D Laplacian */
42: MatCreate(PETSC_COMM_WORLD,&K);
43: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
44: MatSetFromOptions(K);
45: MatSetUp(K);
46: MatGetOwnershipRange(K,&Istart,&Iend);
47: for (II=Istart;II<Iend;II++) {
48: i = II/n; j = II-i*n;
49: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
50: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
51: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
52: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
53: MatSetValue(K,II,II,4.0,INSERT_VALUES);
54: }
55: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
58: /* C is the 1-D Laplacian on horizontal lines */
59: MatCreate(PETSC_COMM_WORLD,&C);
60: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
61: MatSetFromOptions(C);
62: MatSetUp(C);
63: MatGetOwnershipRange(C,&Istart,&Iend);
64: for (II=Istart;II<Iend;II++) {
65: i = II/n; j = II-i*n;
66: if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
67: if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
68: MatSetValue(C,II,II,2.0,INSERT_VALUES);
69: }
70: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
73: /* M is a diagonal matrix */
74: MatCreate(PETSC_COMM_WORLD,&M);
75: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
76: MatSetFromOptions(M);
77: MatSetUp(M);
78: MatGetOwnershipRange(M,&Istart,&Iend);
79: for (II=Istart;II<Iend;II++) {
80: MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
81: }
82: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create a standalone EPS with appropriate settings
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: EPSCreate(PETSC_COMM_WORLD,&eps);
90: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
91: #if defined(PETSC_USE_COMPLEX)
92: EPSSetTarget(eps,0.01*PETSC_i);
93: #endif
94: EPSGetST(eps,&st);
95: STSetType(st,STSINVERT);
96: STGetKSP(st,&ksp);
97: KSPSetType(ksp,KSPBCGS);
98: KSPGetPC(ksp,&pc);
99: PCSetType(pc,PCJACOBI);
100: EPSSetFromOptions(eps);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Create the eigensolver and solve the eigensystem
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106: PEPCreate(PETSC_COMM_WORLD,&pep);
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: PEPLinearGetCompanionForm(pep,&cform);
115: PetscPrintf(PETSC_COMM_WORLD," Linearization with companion form %D",cform);
116: PEPLinearGetExplicitMatrix(pep,&expmat);
117: if (expmat) {
118: PetscPrintf(PETSC_COMM_WORLD," with explicit matrix");
119: }
120: PetscPrintf(PETSC_COMM_WORLD,"\n");
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Display solution and clean up
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
127: PEPDestroy(&pep);
128: EPSDestroy(&eps);
129: MatDestroy(&M);
130: MatDestroy(&C);
131: MatDestroy(&K);
132: SlepcFinalize();
133: return ierr;
134: }