Actual source code: test9.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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[] = "Test the NLEIGS solver with FNCOMBINE.\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],g; /* functions to define the nonlinear operator */
37: NEP nep; /* nonlinear eigensolver context */
38: PetscInt n=100,Istart,Iend,i;
39: PetscReal kappa=1.0,m=1.0;
40: PetscScalar sigma,numer[2],denom[2];
41: PetscBool terse;
44: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
46: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
47: PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL);
48: PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL);
49: sigma = kappa/m;
50: PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%D kappa=%g m=%g\n\n",n,(double)kappa,(double)m);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Build the problem matrices
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: /* initialize matrices */
57: for (i=0;i<NMAT;i++) {
58: MatCreate(PETSC_COMM_WORLD,&A[i]);
59: MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n);
60: MatSetFromOptions(A[i]);
61: MatSetUp(A[i]);
62: }
63: MatGetOwnershipRange(A[0],&Istart,&Iend);
65: /* A0 */
66: for (i=Istart;i<Iend;i++) {
67: MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES);
68: if (i>0) { MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES); }
69: if (i<n-1) { MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES); }
70: }
72: /* A1 */
73: for (i=Istart;i<Iend;i++) {
74: MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES);
75: if (i>0) { MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES); }
76: if (i<n-1) { MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES); }
77: }
79: /* A2 */
80: if (Istart<=n-1 && n-1<Iend) {
81: MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES);
82: }
84: /* assemble matrices */
85: for (i=0;i<NMAT;i++) {
86: MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
87: }
88: for (i=0;i<NMAT;i++) {
89: MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
90: }
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create the problem functions
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /* f1=1 */
97: FNCreate(PETSC_COMM_WORLD,&f[0]);
98: FNSetType(f[0],FNRATIONAL);
99: numer[0] = 1.0;
100: FNRationalSetNumerator(f[0],1,numer);
102: /* f2=-lambda */
103: FNCreate(PETSC_COMM_WORLD,&f[1]);
104: FNSetType(f[1],FNRATIONAL);
105: numer[0] = -1.0; numer[1] = 0.0;
106: FNRationalSetNumerator(f[1],2,numer);
108: /* f3=lambda/(lambda-sigma)=1+sigma/(lambda-sigma) */
109: FNCreate(PETSC_COMM_WORLD,&g);
110: FNSetType(g,FNRATIONAL);
111: numer[0] = sigma;
112: denom[0] = 1.0; denom[1] = -sigma;
113: FNRationalSetNumerator(g,1,numer);
114: FNRationalSetDenominator(g,2,denom);
115: FNCreate(PETSC_COMM_WORLD,&f[2]);
116: FNSetType(f[2],FNCOMBINE);
117: FNCombineSetChildren(f[2],FN_COMBINE_ADD,f[0],g);
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Create the eigensolver and solve the problem
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: NEPCreate(PETSC_COMM_WORLD,&nep);
124: NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN);
125: NEPSetProblemType(nep,NEP_RATIONAL);
126: NEPSetFromOptions(nep);
127: NEPSolve(nep);
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Display solution and clean up
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: /* show detailed info unless -terse option is given by user */
134: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
135: if (terse) {
136: NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
137: } else {
138: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
139: NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
140: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
141: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
142: }
143: NEPDestroy(&nep);
144: for (i=0;i<NMAT;i++) {
145: MatDestroy(&A[i]);
146: FNDestroy(&f[i]);
147: }
148: FNDestroy(&g);
149: SlepcFinalize();
150: return ierr;
151: }
153: /*TEST
155: test:
156: suffix: 1
157: args: -nep_type nleigs -rg_type interval -rg_interval_endpoints 4,700,-.1,.1 -nep_nev 8 -nep_target 5 -terse
158: requires: !single
160: TEST*/