Actual source code: test2.c
slepc-3.20.2 2024-03-15
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[] = "Test NEP interface functions.\n\n";
13: #include <slepcnep.h>
15: int main(int argc,char **argv)
16: {
17: Mat A[3],B; /* problem matrices */
18: FN f[3],g; /* problem functions */
19: NEP nep; /* eigenproblem solver context */
20: DS ds;
21: RG rg;
22: PetscReal tol;
23: PetscScalar coeffs[2],target;
24: PetscInt n=20,i,its,nev,ncv,mpd,Istart,Iend,nterm;
25: PetscBool twoside;
26: NEPWhich which;
27: NEPConvergedReason reason;
28: NEPType type;
29: NEPRefine refine;
30: NEPRefineScheme rscheme;
31: NEPConv conv;
32: NEPStop stop;
33: NEPProblemType ptype;
34: MatStructure mstr;
35: PetscViewerAndFormat *vf;
37: PetscFunctionBeginUser;
38: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
39: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
41: /*
42: Matrices
43: */
44: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
45: PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
46: PetscCall(MatSetFromOptions(A[0]));
47: PetscCall(MatSetUp(A[0]));
48: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
49: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[0],i,i,i+1,INSERT_VALUES));
50: PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
51: PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));
53: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[1]));
54: PetscCall(MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n));
55: PetscCall(MatSetFromOptions(A[1]));
56: PetscCall(MatSetUp(A[1]));
57: PetscCall(MatGetOwnershipRange(A[1],&Istart,&Iend));
58: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[1],i,i,1.0,INSERT_VALUES));
59: PetscCall(MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY));
60: PetscCall(MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY));
62: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[2]));
63: PetscCall(MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n));
64: PetscCall(MatSetFromOptions(A[2]));
65: PetscCall(MatSetUp(A[2]));
66: PetscCall(MatGetOwnershipRange(A[1],&Istart,&Iend));
67: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(A[2],i,i,n/(PetscReal)(i+1),INSERT_VALUES));
68: PetscCall(MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY));
69: PetscCall(MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY));
71: /*
72: Functions: f0=-lambda, f1=1.0, f2=sqrt(lambda)
73: */
74: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
75: PetscCall(FNSetType(f[0],FNRATIONAL));
76: coeffs[0] = -1.0; coeffs[1] = 0.0;
77: PetscCall(FNRationalSetNumerator(f[0],2,coeffs));
79: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
80: PetscCall(FNSetType(f[1],FNRATIONAL));
81: coeffs[0] = 1.0;
82: PetscCall(FNRationalSetNumerator(f[1],1,coeffs));
84: PetscCall(FNCreate(PETSC_COMM_WORLD,&f[2]));
85: PetscCall(FNSetType(f[2],FNSQRT));
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Create eigensolver and test interface functions
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
91: PetscCall(NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN));
92: PetscCall(NEPGetSplitOperatorInfo(nep,&nterm,&mstr));
93: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Nonlinear function with %" PetscInt_FMT " terms, with %s nonzero pattern\n",nterm,MatStructures[mstr]));
94: PetscCall(NEPGetSplitOperatorTerm(nep,0,&B,&g));
95: PetscCall(MatView(B,NULL));
96: PetscCall(FNView(g,NULL));
98: PetscCall(NEPSetType(nep,NEPRII));
99: PetscCall(NEPGetType(nep,&type));
100: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type));
101: PetscCall(NEPGetTwoSided(nep,&twoside));
102: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Two-sided flag = %s\n",twoside?"true":"false"));
104: PetscCall(NEPGetProblemType(nep,&ptype));
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype));
106: PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));
107: PetscCall(NEPGetProblemType(nep,&ptype));
108: PetscCall(PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.\n",(int)ptype));
110: PetscCall(NEPSetRefine(nep,NEP_REFINE_SIMPLE,1,1e-9,2,NEP_REFINE_SCHEME_EXPLICIT));
111: PetscCall(NEPGetRefine(nep,&refine,NULL,&tol,&its,&rscheme));
112: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%" PetscInt_FMT ", scheme=%s\n",NEPRefineTypes[refine],(double)tol,its,NEPRefineSchemes[rscheme]));
114: PetscCall(NEPSetTarget(nep,1.1));
115: PetscCall(NEPGetTarget(nep,&target));
116: PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE));
117: PetscCall(NEPGetWhichEigenpairs(nep,&which));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target)));
120: PetscCall(NEPSetDimensions(nep,1,12,PETSC_DEFAULT));
121: PetscCall(NEPGetDimensions(nep,&nev,&ncv,&mpd));
122: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%" PetscInt_FMT ", ncv=%" PetscInt_FMT ", mpd=%" PetscInt_FMT "\n",nev,ncv,mpd));
124: PetscCall(NEPSetTolerances(nep,1.0e-6,200));
125: PetscCall(NEPGetTolerances(nep,&tol,&its));
126: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.6f, max_its = %" PetscInt_FMT "\n",(double)tol,its));
128: PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_ABS));
129: PetscCall(NEPGetConvergenceTest(nep,&conv));
130: PetscCall(NEPSetStoppingTest(nep,NEP_STOP_BASIC));
131: PetscCall(NEPGetStoppingTest(nep,&stop));
132: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop));
134: PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
135: PetscCall(NEPMonitorSet(nep,(PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))NEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
136: PetscCall(NEPMonitorCancel(nep));
138: PetscCall(NEPGetDS(nep,&ds));
139: PetscCall(DSView(ds,NULL));
140: PetscCall(NEPSetFromOptions(nep));
142: PetscCall(NEPGetRG(nep,&rg));
143: PetscCall(RGView(rg,NULL));
145: PetscCall(NEPSolve(nep));
146: PetscCall(NEPGetConvergedReason(nep,&reason));
147: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason));
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Display solution and clean up
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
153: PetscCall(NEPDestroy(&nep));
154: PetscCall(MatDestroy(&A[0]));
155: PetscCall(MatDestroy(&A[1]));
156: PetscCall(MatDestroy(&A[2]));
157: PetscCall(FNDestroy(&f[0]));
158: PetscCall(FNDestroy(&f[1]));
159: PetscCall(FNDestroy(&f[2]));
160: PetscCall(SlepcFinalize());
161: return 0;
162: }
164: /*TEST
166: test:
167: suffix: 1
169: TEST*/