Actual source code: test2.c
slepc-3.13.0 2020-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: PetscErrorCode ierr;
36: PetscViewerAndFormat *vf;
37: const char *str=NULL;
39: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
40: PetscPrintf(PETSC_COMM_WORLD,"\nDiagonal Nonlinear Eigenproblem, n=%D\n\n",n);
42: /*
43: Matrices
44: */
45: MatCreate(PETSC_COMM_WORLD,&A[0]);
46: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
47: MatSetFromOptions(A[0]);
48: MatSetUp(A[0]);
49: MatGetOwnershipRange(A[0],&Istart,&Iend);
50: for (i=Istart;i<Iend;i++) {
51: MatSetValue(A[0],i,i,i+1,INSERT_VALUES);
52: }
53: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
56: MatCreate(PETSC_COMM_WORLD,&A[1]);
57: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(A[1]);
59: MatSetUp(A[1]);
60: MatGetOwnershipRange(A[1],&Istart,&Iend);
61: for (i=Istart;i<Iend;i++) {
62: MatSetValue(A[1],i,i,1.0,INSERT_VALUES);
63: }
64: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
67: MatCreate(PETSC_COMM_WORLD,&A[2]);
68: MatSetSizes(A[2],PETSC_DECIDE,PETSC_DECIDE,n,n);
69: MatSetFromOptions(A[2]);
70: MatSetUp(A[2]);
71: MatGetOwnershipRange(A[1],&Istart,&Iend);
72: for (i=Istart;i<Iend;i++) {
73: MatSetValue(A[2],i,i,n/(PetscReal)(i+1),INSERT_VALUES);
74: }
75: MatAssemblyBegin(A[2],MAT_FINAL_ASSEMBLY);
76: MatAssemblyEnd(A[2],MAT_FINAL_ASSEMBLY);
78: /*
79: Functions: f0=-lambda, f1=1.0, f2=sqrt(lambda)
80: */
81: FNCreate(PETSC_COMM_WORLD,&f[0]);
82: FNSetType(f[0],FNRATIONAL);
83: coeffs[0] = -1.0; coeffs[1] = 0.0;
84: FNRationalSetNumerator(f[0],2,coeffs);
86: FNCreate(PETSC_COMM_WORLD,&f[1]);
87: FNSetType(f[1],FNRATIONAL);
88: coeffs[0] = 1.0;
89: FNRationalSetNumerator(f[1],1,coeffs);
91: FNCreate(PETSC_COMM_WORLD,&f[2]);
92: FNSetType(f[2],FNSQRT);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create eigensolver and test interface functions
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: NEPCreate(PETSC_COMM_WORLD,&nep);
98: NEPSetSplitOperator(nep,3,A,f,SAME_NONZERO_PATTERN);
99: NEPGetSplitOperatorInfo(nep,&nterm,&mstr);
100: switch (mstr) {
101: case DIFFERENT_NONZERO_PATTERN: str = "different"; break;
102: case SUBSET_NONZERO_PATTERN: str = "subset"; break;
103: case SAME_NONZERO_PATTERN: str = "same"; break;
104: }
105: PetscPrintf(PETSC_COMM_WORLD," Nonlinear function with %d terms, with %s nonzero pattern\n",(int)nterm,str);
106: NEPGetSplitOperatorTerm(nep,0,&B,&g);
107: MatView(B,NULL);
108: FNView(g,NULL);
110: NEPSetType(nep,NEPRII);
111: NEPGetType(nep,&type);
112: PetscPrintf(PETSC_COMM_WORLD," Type set to %s\n",type);
113: NEPGetTwoSided(nep,&twoside);
114: PetscPrintf(PETSC_COMM_WORLD," Two-sided flag = %s\n",twoside?"true":"false");
116: NEPGetProblemType(nep,&ptype);
117: PetscPrintf(PETSC_COMM_WORLD," Problem type before changing = %d",(int)ptype);
118: NEPSetProblemType(nep,NEP_RATIONAL);
119: NEPGetProblemType(nep,&ptype);
120: PetscPrintf(PETSC_COMM_WORLD," ... changed to %d.\n",(int)ptype);
122: NEPSetRefine(nep,NEP_REFINE_SIMPLE,1,1e-9,2,NEP_REFINE_SCHEME_EXPLICIT);
123: NEPGetRefine(nep,&refine,NULL,&tol,&its,&rscheme);
124: PetscPrintf(PETSC_COMM_WORLD," Refinement: %s, tol=%g, its=%D, scheme=%s\n",NEPRefineTypes[refine],(double)tol,its,NEPRefineSchemes[rscheme]);
126: NEPSetTarget(nep,1.1);
127: NEPGetTarget(nep,&target);
128: NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
129: NEPGetWhichEigenpairs(nep,&which);
130: PetscPrintf(PETSC_COMM_WORLD," Which = %d, target = %g\n",(int)which,(double)PetscRealPart(target));
132: NEPSetDimensions(nep,1,12,PETSC_DEFAULT);
133: NEPGetDimensions(nep,&nev,&ncv,&mpd);
134: PetscPrintf(PETSC_COMM_WORLD," Dimensions: nev=%D, ncv=%D, mpd=%D\n",nev,ncv,mpd);
136: NEPSetTolerances(nep,1.0e-6,200);
137: NEPGetTolerances(nep,&tol,&its);
138: PetscPrintf(PETSC_COMM_WORLD," Tolerance = %.6f, max_its = %D\n",(double)tol,its);
140: NEPSetConvergenceTest(nep,NEP_CONV_ABS);
141: NEPGetConvergenceTest(nep,&conv);
142: NEPSetStoppingTest(nep,NEP_STOP_BASIC);
143: NEPGetStoppingTest(nep,&stop);
144: PetscPrintf(PETSC_COMM_WORLD," Convergence test = %d, stopping test = %d\n",(int)conv,(int)stop);
146: PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
147: NEPMonitorSet(nep,(PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))NEPMonitorFirst,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
148: NEPMonitorCancel(nep);
150: NEPGetDS(nep,&ds);
151: DSView(ds,NULL);
152: NEPSetFromOptions(nep);
154: NEPGetRG(nep,&rg);
155: RGView(rg,NULL);
157: NEPSolve(nep);
158: NEPGetConvergedReason(nep,&reason);
159: PetscPrintf(PETSC_COMM_WORLD," Finished - converged reason = %d\n",(int)reason);
161: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162: Display solution and clean up
163: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164: NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
165: NEPDestroy(&nep);
166: MatDestroy(&A[0]);
167: MatDestroy(&A[1]);
168: MatDestroy(&A[2]);
169: FNDestroy(&f[0]);
170: FNDestroy(&f[1]);
171: FNDestroy(&f[2]);
172: SlepcFinalize();
173: return ierr;
174: }
176: /*TEST
178: test:
179: suffix: 1
180: requires: !single
182: TEST*/