Actual source code: test12.c
slepc-3.12.2 2020-01-13
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: */
11: static char help[] = "Test some NLEIGS interface functions.\n\n"
12: "Based on ex27.c. The command line options are:\n"
13: " -n <n>, where <n> = matrix dimension.\n";
15: /*
16: Solve T(lambda)x=0 using NLEIGS solver
17: with T(lambda) = -D+sqrt(lambda)*I
18: where D is the Laplacian operator in 1 dimension
19: and with the interpolation interval [.01,16]
20: */
22: #include <slepcnep.h>
24: /*
25: User-defined routines
26: */
27: PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);
29: int main(int argc,char **argv)
30: {
31: NEP nep; /* nonlinear eigensolver context */
32: Mat A[2];
33: PetscInt n=100,Istart,Iend,i,ns,nsin;
35: PetscBool terse,fb;
36: RG rg;
37: FN f[2];
38: PetscScalar coeffs,shifts[]={1.06,1.1,1.12,1.15},*rkshifts,val;
39: PetscErrorCode (*fsing)(NEP,PetscInt*,PetscScalar*,void*);
41: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
43: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%D\n\n",n);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Create nonlinear eigensolver and set some options
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: NEPCreate(PETSC_COMM_WORLD,&nep);
50: NEPSetType(nep,NEPNLEIGS);
51: NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL);
52: NEPGetRG(nep,&rg);
53: RGSetType(rg,RGINTERVAL);
54: #if defined(PETSC_USE_COMPLEX)
55: RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001);
56: #else
57: RGIntervalSetEndpoints(rg,0.01,16.0,0,0);
58: #endif
59: NEPSetTarget(nep,1.1);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Define the nonlinear problem in split form
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: /* Create matrices */
66: MatCreate(PETSC_COMM_WORLD,&A[0]);
67: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
68: MatSetFromOptions(A[0]);
69: MatSetUp(A[0]);
70: MatGetOwnershipRange(A[0],&Istart,&Iend);
71: for (i=Istart;i<Iend;i++) {
72: if (i>0) { MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES); }
73: if (i<n-1) { MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES); }
74: MatSetValue(A[0],i,i,-2.0,INSERT_VALUES);
75: }
76: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
79: MatCreate(PETSC_COMM_WORLD,&A[1]);
80: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
81: MatSetFromOptions(A[1]);
82: MatSetUp(A[1]);
83: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
85: MatShift(A[1],1.0);
87: /* Define functions */
88: FNCreate(PETSC_COMM_WORLD,&f[0]);
89: FNSetType(f[0],FNRATIONAL);
90: coeffs = 1.0;
91: FNRationalSetNumerator(f[0],1,&coeffs);
92: FNCreate(PETSC_COMM_WORLD,&f[1]);
93: FNSetType(f[1],FNSQRT);
94: NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Set some options
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: NEPNLEIGSSetFullBasis(nep,PETSC_FALSE);
101: NEPNLEIGSSetRKShifts(nep,4,shifts);
102: NEPSetFromOptions(nep);
104: NEPNLEIGSGetFullBasis(nep,&fb);
105: PetscPrintf(PETSC_COMM_WORLD," Using full basis = %s\n",fb?"true":"false");
106: NEPNLEIGSGetRKShifts(nep,&ns,&rkshifts);
107: if (ns) {
108: PetscPrintf(PETSC_COMM_WORLD," Using %d RK shifts =",ns);
109: for (i=0;i<ns;i++) {
110: PetscPrintf(PETSC_COMM_WORLD," %g",(double)PetscRealPart(rkshifts[i]));
111: }
112: PetscPrintf(PETSC_COMM_WORLD,"\n");
113: PetscFree(rkshifts);
114: }
115: NEPNLEIGSGetSingularitiesFunction(nep,&fsing,NULL);
116: nsin = 1;
117: (*fsing)(nep,&nsin,&val,NULL);
118: PetscPrintf(PETSC_COMM_WORLD," First returned singularity = %g\n",(double)PetscRealPart(val));
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Solve the eigensystem
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: NEPSolve(nep);
125: /* show detailed info unless -terse option is given by user */
126: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
127: if (terse) {
128: NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
129: } else {
130: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
131: NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
132: NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
133: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
134: }
135: NEPDestroy(&nep);
136: MatDestroy(&A[0]);
137: MatDestroy(&A[1]);
138: FNDestroy(&f[0]);
139: FNDestroy(&f[1]);
140: SlepcFinalize();
141: return ierr;
142: }
144: /* ------------------------------------------------------------------- */
145: /*
146: ComputeSingularities - Computes maxnp points (at most) in the complex plane where
147: the function T(.) is not analytic.
149: In this case, we discretize the singularity region (-inf,0)~(-10e+6,-10e-6)
150: */
151: PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
152: {
153: PetscReal h;
154: PetscInt i;
157: h = 11.0/(*maxnp-1);
158: xi[0] = -1e-5; xi[*maxnp-1] = -1e+6;
159: for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-5+h*i);
160: return(0);
161: }
163: /*TEST
165: test:
166: suffix: 1
167: args: -nep_nev 3 -nep_nleigs_interpolation_degree 90 -terse -nep_view
168: requires: double
169: filter: grep -v tolerance
171: TEST*/