Actual source code: test9.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: */
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;
43: PetscFunctionBeginUser;
44: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
46: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
47: PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
48: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mass",&m,NULL));
49: sigma = kappa/m;
50: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Loaded vibrating string, n=%" PetscInt_FMT " 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: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
59: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
60: PetscCall(MatSetFromOptions(A[i]));
61: PetscCall(MatSetUp(A[i]));
62: }
63: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
65: /* A0 */
66: for (i=Istart;i<Iend;i++) {
67: PetscCall(MatSetValue(A[0],i,i,(i==n-1)?1.0*n:2.0*n,INSERT_VALUES));
68: if (i>0) PetscCall(MatSetValue(A[0],i,i-1,-1.0*n,INSERT_VALUES));
69: if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,-1.0*n,INSERT_VALUES));
70: }
72: /* A1 */
73: for (i=Istart;i<Iend;i++) {
74: PetscCall(MatSetValue(A[1],i,i,(i==n-1)?2.0/(6.0*n):4.0/(6.0*n),INSERT_VALUES));
75: if (i>0) PetscCall(MatSetValue(A[1],i,i-1,1.0/(6.0*n),INSERT_VALUES));
76: if (i<n-1) PetscCall(MatSetValue(A[1],i,i+1,1.0/(6.0*n),INSERT_VALUES));
77: }
79: /* A2 */
80: if (Istart<=n-1 && n-1<Iend) PetscCall(MatSetValue(A[2],n-1,n-1,kappa,INSERT_VALUES));
82: /* assemble matrices */
83: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
84: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Create the problem functions
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: /* f1=1 */
91: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
92: PetscCall(FNSetType(f[0],FNRATIONAL));
93: numer[0] = 1.0;
94: PetscCall(FNRationalSetNumerator(f[0],1,numer));
96: /* f2=-lambda */
97: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
98: PetscCall(FNSetType(f[1],FNRATIONAL));
99: numer[0] = -1.0; numer[1] = 0.0;
100: PetscCall(FNRationalSetNumerator(f[1],2,numer));
102: /* f3=lambda/(lambda-sigma)=1+sigma/(lambda-sigma) */
103: PetscCall(FNCreate(PETSC_COMM_WORLD,&g));
104: PetscCall(FNSetType(g,FNRATIONAL));
105: numer[0] = sigma;
106: denom[0] = 1.0; denom[1] = -sigma;
107: PetscCall(FNRationalSetNumerator(g,1,numer));
108: PetscCall(FNRationalSetDenominator(g,2,denom));
109: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
110: PetscCall(FNSetType(f[2],FNCOMBINE));
111: PetscCall(FNCombineSetChildren(f[2],FN_COMBINE_ADD,f[0],g));
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Create the eigensolver and solve the problem
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
118: PetscCall(NEPSetSplitOperator(nep,3,A,f,SUBSET_NONZERO_PATTERN));
119: PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
120: PetscCall(NEPSetFromOptions(nep));
121: PetscCall(NEPSolve(nep));
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Display solution and clean up
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: /* show detailed info unless -terse option is given by user */
128: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
129: if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
130: else {
131: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
132: PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
133: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
134: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
135: }
136: PetscCall(NEPDestroy(&nep));
137: for (i=0;i<NMAT;i++) {
138: PetscCall(MatDestroy(&A[i]));
139: PetscCall(FNDestroy(&f[i]));
140: }
141: PetscCall(FNDestroy(&g));
142: PetscCall(SlepcFinalize());
143: return 0;
144: }
146: /*TEST
148: test:
149: suffix: 1
150: args: -nep_type nleigs -rg_type interval -rg_interval_endpoints 4,700,-.1,.1 -nep_nev 8 -nep_target 5 -terse
151: filter: sed -e "s/[+-]0\.0*i//g"
152: requires: !single
154: TEST*/