Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test DSNEP.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 2 : int main(int argc,char **argv)
16 : {
17 2 : DS ds;
18 2 : FN f1,f2,f3,funs[3],qfun;
19 2 : SlepcSC sc;
20 2 : PetscScalar *Id,*A,*B,*wr,*wi,*X,*W,coeffs[2],auxr,alpha;
21 2 : PetscReal tol,tau=0.001,radius=10,h,a=20,xi,re,im,nrm,aux;
22 2 : PetscInt i,j,ii,jj,k,n=10,ld,nev,nfun;
23 2 : PetscViewer viewer;
24 2 : PetscBool verbose;
25 2 : RG rg;
26 2 : DSMatType mat[3]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2};
27 : #if !defined(PETSC_USE_COMPLEX)
28 : PetscScalar auxi;
29 : #endif
30 :
31 2 : PetscFunctionBeginUser;
32 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
33 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
34 2 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
35 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %" PetscInt_FMT ", tau=%g.\n",n,(double)tau));
36 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
37 2 : PetscCall(PetscOptionsGetReal(NULL,NULL,"-radius",&radius,NULL));
38 :
39 : /* Create DS object */
40 2 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
41 2 : PetscCall(DSSetType(ds,DSNEP));
42 2 : tol = 1000*n*PETSC_MACHINE_EPSILON;
43 2 : PetscCall(DSNEPSetRefine(ds,tol,PETSC_DECIDE));
44 2 : PetscCall(DSSetFromOptions(ds));
45 :
46 : /* Set functions (prior to DSAllocate) */
47 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
48 2 : PetscCall(FNSetType(f1,FNRATIONAL));
49 2 : coeffs[0] = -1.0; coeffs[1] = 0.0;
50 2 : PetscCall(FNRationalSetNumerator(f1,2,coeffs));
51 :
52 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
53 2 : PetscCall(FNSetType(f2,FNRATIONAL));
54 2 : coeffs[0] = 1.0;
55 2 : PetscCall(FNRationalSetNumerator(f2,1,coeffs));
56 :
57 2 : PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
58 2 : PetscCall(FNSetType(f3,FNEXP));
59 2 : PetscCall(FNSetScale(f3,-tau,1.0));
60 :
61 2 : funs[0] = f1;
62 2 : funs[1] = f2;
63 2 : funs[2] = f3;
64 2 : PetscCall(DSNEPSetFN(ds,3,funs));
65 :
66 : /* Set dimensions */
67 2 : ld = n+2; /* test leading dimension larger than n */
68 2 : PetscCall(DSAllocate(ds,ld));
69 2 : PetscCall(DSSetDimensions(ds,n,0,0));
70 :
71 : /* Set region (used only in method=1) */
72 2 : PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
73 2 : PetscCall(RGSetType(rg,RGELLIPSE));
74 2 : PetscCall(RGEllipseSetParameters(rg,0.0,radius,1.0));
75 2 : PetscCall(DSNEPSetRG(ds,rg));
76 2 : PetscCall(RGDestroy(&rg));
77 :
78 : /* Set up viewer */
79 2 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
80 2 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
81 2 : PetscCall(DSView(ds,viewer));
82 2 : PetscCall(PetscViewerPopFormat(viewer));
83 2 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
84 :
85 : /* Show info about functions */
86 2 : PetscCall(DSNEPGetNumFN(ds,&nfun));
87 8 : for (i=0;i<nfun;i++) {
88 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Function %" PetscInt_FMT ":\n",i));
89 6 : PetscCall(DSNEPGetFN(ds,i,&qfun));
90 6 : PetscCall(FNView(qfun,NULL));
91 : }
92 :
93 : /* Fill matrices */
94 2 : PetscCall(DSGetArray(ds,DS_MAT_E0,&Id));
95 22 : for (i=0;i<n;i++) Id[i+i*ld]=1.0;
96 2 : PetscCall(DSRestoreArray(ds,DS_MAT_E0,&Id));
97 2 : h = PETSC_PI/(PetscReal)(n+1);
98 2 : PetscCall(DSGetArray(ds,DS_MAT_E1,&A));
99 22 : for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
100 20 : for (i=1;i<n;i++) {
101 18 : A[i+(i-1)*ld]=1.0/(h*h);
102 18 : A[(i-1)+i*ld]=1.0/(h*h);
103 : }
104 2 : PetscCall(DSRestoreArray(ds,DS_MAT_E1,&A));
105 2 : PetscCall(DSGetArray(ds,DS_MAT_E2,&B));
106 22 : for (i=0;i<n;i++) {
107 20 : xi = (i+1)*h;
108 20 : B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
109 : }
110 2 : PetscCall(DSRestoreArray(ds,DS_MAT_E2,&B));
111 :
112 2 : if (verbose) {
113 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
114 0 : PetscCall(DSView(ds,viewer));
115 : }
116 :
117 : /* Solve */
118 2 : PetscCall(PetscCalloc2(n,&wr,n,&wi));
119 2 : PetscCall(DSGetSlepcSC(ds,&sc));
120 2 : sc->comparison = SlepcCompareLargestMagnitude;
121 2 : sc->comparisonctx = NULL;
122 2 : sc->map = NULL;
123 2 : sc->mapobj = NULL;
124 2 : PetscCall(DSSolve(ds,wr,wi));
125 2 : PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));
126 :
127 2 : if (verbose) {
128 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
129 0 : PetscCall(DSView(ds,viewer));
130 : }
131 2 : PetscCall(DSGetDimensions(ds,NULL,NULL,NULL,&nev));
132 :
133 : /* Print computed eigenvalues */
134 2 : PetscCall(PetscMalloc1(ld*ld,&W));
135 2 : PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL));
136 2 : PetscCall(DSGetArray(ds,DS_MAT_X,&X));
137 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
138 6 : for (i=0;i<nev;i++) {
139 : #if defined(PETSC_USE_COMPLEX)
140 4 : re = PetscRealPart(wr[i]);
141 4 : im = PetscImaginaryPart(wr[i]);
142 : #else
143 : re = wr[i];
144 : im = wi[i];
145 : #endif
146 : /* Residual */
147 4 : PetscCall(PetscArrayzero(W,ld*ld));
148 16 : for (k=0;k<nfun;k++) {
149 12 : PetscCall(FNEvaluateFunction(funs[k],wr[i],&alpha));
150 12 : PetscCall(DSGetArray(ds,mat[k],&A));
151 1332 : for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
152 12 : PetscCall(DSRestoreArray(ds,mat[k],&A));
153 : }
154 : nrm = 0.0;
155 44 : for (k=0;k<n;k++) {
156 : auxr = 0.0;
157 : #if !defined(PETSC_USE_COMPLEX)
158 : auxi = 0.0;
159 : #endif
160 440 : for (j=0;j<n;j++) {
161 400 : auxr += W[k+j*ld]*X[i*ld+j];
162 : #if !defined(PETSC_USE_COMPLEX)
163 : if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
164 : #endif
165 : }
166 40 : aux = SlepcAbsEigenvalue(auxr,auxi);
167 40 : nrm += aux*aux;
168 : }
169 4 : nrm = PetscSqrtReal(nrm);
170 4 : if (nrm/SlepcAbsEigenvalue(wr[i],wi[i])>tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %" PetscInt_FMT "-th computed eigenpair %g\n",i,(double)nrm));
171 4 : if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
172 4 : else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
173 : }
174 2 : PetscCall(PetscFree(W));
175 2 : PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
176 2 : PetscCall(DSRestoreArray(ds,DS_MAT_W,&W));
177 2 : PetscCall(PetscFree2(wr,wi));
178 2 : PetscCall(FNDestroy(&f1));
179 2 : PetscCall(FNDestroy(&f2));
180 2 : PetscCall(FNDestroy(&f3));
181 2 : PetscCall(DSDestroy(&ds));
182 2 : PetscCall(SlepcFinalize());
183 : return 0;
184 : }
185 :
186 : /*TEST
187 :
188 : testset:
189 : test:
190 : filter: grep -v "solving the problem"
191 : suffix: 1
192 : test:
193 : suffix: 2
194 : args: -ds_method 1 -radius 10 -ds_nep_refine_its 1
195 : filter: grep -v "solving the problem" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/" | sed -e "s/tolerance [0-9]\.[0-9]*e[+-]\([0-9]*\)/tolerance removed/" | sed -e "s/tolerance [0-9]\.\([0-9]*\)/tolerance removed/"
196 : requires: complex
197 :
198 : TEST*/
|