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