Actual source code: ex42.c
slepc-3.20.1 2023-11-27
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: */
10: /*
11: This example implements one of the problems found at
12: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
13: The University of Manchester.
14: The details of the collection can be found at:
15: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
16: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
18: The loaded_string problem is a rational eigenvalue problem for the
19: finite element model of a loaded vibrating string.
20: */
22: static char help[] = "Illustrates computation of left eigenvectors and resolvent.\n\n"
23: "This is based on loaded_string from the NLEVP collection.\n"
24: "The command line options are:\n"
25: " -n <n>, dimension of the matrices.\n"
26: " -kappa <kappa>, stiffness of elastic spring.\n"
27: " -mass <m>, mass of the attached load.\n\n";
29: #include <slepcnep.h>
31: #define NMAT 3
33: int main(int argc,char **argv)
34: {
35: Mat A[NMAT]; /* problem matrices */
36: FN f[NMAT]; /* functions to define the nonlinear operator */
37: NEP nep; /* nonlinear eigensolver context */
38: RG rg;
39: Vec v,r,z,w;
40: PetscInt n=100,Istart,Iend,i,nconv;
41: PetscReal kappa=1.0,m=1.0,nrm,tol;
42: PetscScalar lambda,sigma,numer[2],denom[2],omega1,omega2;
44: PetscFunctionBeginUser;
45: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
47: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
48: PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
49: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
50: sigma = kappa/m;
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%" PetscInt_FMT " kappa=%g m=%g\n\n",n,(double)kappa,(double)m));
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Build the problem matrices
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: /* initialize matrices */
58: for (i=0;i<NMAT;i++) {
59: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
60: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
61: PetscCall(MatSetFromOptions(A[i]));
62: PetscCall(MatSetUp(A[i]));
63: }
64: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
66: /* A0 */
67: for (i=Istart;i<Iend;i++) {
68: PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
69: if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
70: if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
71: }
73: /* A1 */
74: for (i=Istart;i<Iend;i++) {
75: PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
76: if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
77: if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
78: }
80: /* A2 */
81: if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));
83: /* assemble matrices */
84: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
85: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Create the problem functions
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: /* f1=1 */
92: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
93: PetscCall(FNSetType(f[0],FNRATIONAL));
94: numer[0] = 1.0;
95: PetscCall(FNRationalSetNumerator(f[0],1,numer));
97: /* f2=-lambda */
98: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
99: PetscCall(FNSetType(f[1],FNRATIONAL));
100: numer[0] = -1.0; numer[1] = 0.0;
101: PetscCall(FNRationalSetNumerator(f[1],2,numer));
103: /* f3=lambda/(lambda-sigma) */
104: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
105: PetscCall(FNSetType(f[2],FNRATIONAL));
106: numer[0] = 1.0; numer[1] = 0.0;
107: denom[0] = 1.0; denom[1] = -sigma;
108: PetscCall(FNRationalSetNumerator(f[2],2,numer));
109: PetscCall(FNRationalSetDenominator(f[2],2,denom));
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Create the eigensolver and solve the problem
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
116: PetscCall(NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN));
117: PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
118: PetscCall(NEPSetDimensions(nep,8,PETSC_DEFAULT,PETSC_DEFAULT));
120: /* set two-sided NLEIGS solver */
121: PetscCall(NEPSetType(nep,NEPNLEIGS));
122: PetscCall(NEPNLEIGSSetFullBasis(nep,PETSC_TRUE));
123: PetscCall(NEPSetTwoSided(nep,PETSC_TRUE));
124: PetscCall(NEPGetRG(nep,&rg));
125: PetscCall(RGSetType(rg,RGINTERVAL));
126: #if defined(PETSC_USE_COMPLEX)
127: PetscCall(RGIntervalSetEndpoints(rg,4.0,700.0,-0.001,0.001));
128: #else
129: PetscCall(RGIntervalSetEndpoints(rg,4.0,700.0,0,0));
130: #endif
131: PetscCall(NEPSetTarget(nep,5.0));
133: PetscCall(NEPSetFromOptions(nep));
134: PetscCall(NEPSolve(nep));
136: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137: Check left residual
138: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139: PetscCall(MatCreateVecs(A[0],&v,&r));
140: PetscCall(VecDuplicate(v,&w));
141: PetscCall(VecDuplicate(v,&z));
142: PetscCall(NEPGetConverged(nep,&nconv));
143: PetscCall(NEPGetTolerances(nep,&tol,NULL));
144: for (i=0;i<nconv;i++) {
145: PetscCall(NEPGetEigenpair(nep,i,&lambda,NULL,NULL,NULL));
146: PetscCall(NEPGetLeftEigenvector(nep,i,v,NULL));
147: PetscCall(NEPApplyAdjoint(nep,lambda,v,w,r,NULL,NULL));
148: PetscCall(VecNorm(r,NORM_2,&nrm));
149: if (nrm>tol*PetscAbsScalar(lambda)) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Left residual i=%" PetscInt_FMT " is above tolerance --> %g\n",i,(double)(nrm/PetscAbsScalar(lambda))));
150: }
152: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153: Operate with resolvent
154: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: omega1 = 20.0;
156: omega2 = 150.0;
157: PetscCall(VecSet(v,0.0));
158: PetscCall(VecSetValue(v,0,-1.0,INSERT_VALUES));
159: PetscCall(VecSetValue(v,1,3.0,INSERT_VALUES));
160: PetscCall(VecAssemblyBegin(v));
161: PetscCall(VecAssemblyEnd(v));
162: PetscCall(NEPApplyResolvent(nep,NULL,omega1,v,r));
163: PetscCall(VecNorm(r,NORM_2,&nrm));
164: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega1),(double)nrm));
165: PetscCall(NEPApplyResolvent(nep,NULL,omega2,v,r));
166: PetscCall(VecNorm(r,NORM_2,&nrm));
167: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega2),(double)nrm));
168: PetscCall(VecSet(v,1.0));
169: PetscCall(NEPApplyResolvent(nep,NULL,omega1,v,r));
170: PetscCall(VecNorm(r,NORM_2,&nrm));
171: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega1),(double)nrm));
172: PetscCall(NEPApplyResolvent(nep,NULL,omega2,v,r));
173: PetscCall(VecNorm(r,NORM_2,&nrm));
174: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"resolvent, omega=%g: norm of computed vector=%g\n",(double)PetscRealPart(omega2),(double)nrm));
176: /* clean up */
177: PetscCall(NEPDestroy(&nep));
178: for (i=0;i<NMAT;i++) {
179: PetscCall(MatDestroy(&A[i]));
180: PetscCall(FNDestroy(&f[i]));
181: }
182: PetscCall(VecDestroy(&v));
183: PetscCall(VecDestroy(&r));
184: PetscCall(VecDestroy(&w));
185: PetscCall(VecDestroy(&z));
186: PetscCall(SlepcFinalize());
187: return 0;
188: }
190: /*TEST
192: test:
193: suffix: 1
194: requires: !single
196: TEST*/