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 RG interface functions.\n\n";
12 :
13 : #include <slepcrg.h>
14 :
15 : #define NPOINTS 10
16 : #define NVERTEX 7
17 :
18 2 : int main(int argc,char **argv)
19 : {
20 2 : RG rg;
21 2 : PetscInt i,inside,nv;
22 2 : PetscBool triv;
23 2 : PetscReal re,im,a,b,c,d;
24 2 : PetscScalar ar,ai,cr[NPOINTS],ci[NPOINTS],vr[NVERTEX],vi[NVERTEX],*pr,*pi;
25 :
26 2 : PetscFunctionBeginUser;
27 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
28 2 : PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
29 :
30 : /* ellipse */
31 2 : PetscCall(RGSetType(rg,RGELLIPSE));
32 2 : PetscCall(RGIsTrivial(rg,&triv));
33 2 : PetscCheck(triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be trivial before setting parameters");
34 2 : PetscCall(RGEllipseSetParameters(rg,1.1,2,0.1));
35 2 : PetscCall(RGSetFromOptions(rg));
36 2 : PetscCall(RGIsTrivial(rg,&triv));
37 2 : PetscCheck(!triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be non-trivial after setting parameters");
38 2 : PetscCall(RGView(rg,NULL));
39 2 : PetscCall(RGViewFromOptions(rg,NULL,"-rg_ellipse_view"));
40 2 : re = 0.1; im = 0.3;
41 : #if defined(PETSC_USE_COMPLEX)
42 : ar = PetscCMPLX(re,im);
43 : #else
44 2 : ar = re; ai = im;
45 : #endif
46 2 : PetscCall(RGCheckInside(rg,1,&ar,&ai,&inside));
47 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point (%g,%g) is %s the region\n",(double)re,(double)im,(inside>=0)?"inside":"outside"));
48 :
49 2 : PetscCall(RGComputeBoundingBox(rg,&a,&b,&c,&d));
50 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The bounding box is [%g,%g]x[%g,%g]\n",(double)a,(double)b,(double)c,(double)d));
51 :
52 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Contour points: "));
53 2 : PetscCall(RGComputeContour(rg,NPOINTS,cr,ci));
54 22 : for (i=0;i<NPOINTS;i++) {
55 : #if defined(PETSC_USE_COMPLEX)
56 : re = PetscRealPart(cr[i]);
57 : im = PetscImaginaryPart(cr[i]);
58 : #else
59 20 : re = cr[i];
60 20 : im = ci[i];
61 : #endif
62 20 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"(%.3g,%.3g) ",(double)re,(double)im));
63 : }
64 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
65 :
66 : /* interval */
67 2 : PetscCall(RGSetType(rg,RGINTERVAL));
68 2 : PetscCall(RGIsTrivial(rg,&triv));
69 2 : PetscCheck(triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be trivial before setting parameters");
70 2 : PetscCall(RGIntervalSetEndpoints(rg,-1,1,-0.1,0.1));
71 2 : PetscCall(RGSetFromOptions(rg));
72 2 : PetscCall(RGIsTrivial(rg,&triv));
73 2 : PetscCheck(!triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be non-trivial after setting parameters");
74 2 : PetscCall(RGView(rg,NULL));
75 2 : PetscCall(RGViewFromOptions(rg,NULL,"-rg_interval_view"));
76 2 : re = 0.2; im = 0;
77 : #if defined(PETSC_USE_COMPLEX)
78 : ar = PetscCMPLX(re,im);
79 : #else
80 2 : ar = re; ai = im;
81 : #endif
82 2 : PetscCall(RGCheckInside(rg,1,&ar,&ai,&inside));
83 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point (%g,%g) is %s the region\n",(double)re,(double)im,(inside>=0)?"inside":"outside"));
84 :
85 2 : PetscCall(RGComputeBoundingBox(rg,&a,&b,&c,&d));
86 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The bounding box is [%g,%g]x[%g,%g]\n",(double)a,(double)b,(double)c,(double)d));
87 :
88 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Contour points: "));
89 2 : PetscCall(RGComputeContour(rg,NPOINTS,cr,ci));
90 22 : for (i=0;i<NPOINTS;i++) {
91 : #if defined(PETSC_USE_COMPLEX)
92 : re = PetscRealPart(cr[i]);
93 : im = PetscImaginaryPart(cr[i]);
94 : #else
95 20 : re = cr[i];
96 20 : im = ci[i];
97 : #endif
98 20 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"(%.3g,%.3g) ",(double)re,(double)im));
99 : }
100 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
101 :
102 : /* polygon */
103 : #if defined(PETSC_USE_COMPLEX)
104 : vr[0] = PetscCMPLX(0.0,2.0);
105 : vr[1] = PetscCMPLX(1.0,4.0);
106 : vr[2] = PetscCMPLX(2.0,5.0);
107 : vr[3] = PetscCMPLX(4.0,3.0);
108 : vr[4] = PetscCMPLX(5.0,4.0);
109 : vr[5] = PetscCMPLX(6.0,1.0);
110 : vr[6] = PetscCMPLX(2.0,0.0);
111 : #else
112 2 : vr[0] = 0.0; vi[0] = 1.0;
113 2 : vr[1] = 0.0; vi[1] = -1.0;
114 2 : vr[2] = 0.6; vi[2] = -0.8;
115 2 : vr[3] = 1.0; vi[3] = -1.0;
116 2 : vr[4] = 2.0; vi[4] = 0.0;
117 2 : vr[5] = 1.0; vi[5] = 1.0;
118 2 : vr[6] = 0.6; vi[6] = 0.8;
119 : #endif
120 2 : PetscCall(RGSetType(rg,RGPOLYGON));
121 2 : PetscCall(RGIsTrivial(rg,&triv));
122 2 : PetscCheck(triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be trivial before setting parameters");
123 2 : PetscCall(RGPolygonSetVertices(rg,NVERTEX,vr,vi));
124 2 : PetscCall(RGSetFromOptions(rg));
125 2 : PetscCall(RGIsTrivial(rg,&triv));
126 2 : PetscCheck(!triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be non-trivial after setting parameters");
127 2 : PetscCall(RGView(rg,NULL));
128 2 : PetscCall(RGViewFromOptions(rg,NULL,"-rg_polygon_view"));
129 2 : re = 5; im = 0.9;
130 : #if defined(PETSC_USE_COMPLEX)
131 : ar = PetscCMPLX(re,im);
132 : #else
133 2 : ar = re; ai = im;
134 : #endif
135 2 : PetscCall(RGCheckInside(rg,1,&ar,&ai,&inside));
136 4 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point (%g,%g) is %s the region\n",(double)re,(double)im,(inside>=0)?"inside":"outside"));
137 :
138 2 : PetscCall(RGComputeBoundingBox(rg,&a,&b,&c,&d));
139 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The bounding box is [%g,%g]x[%g,%g]\n",(double)a,(double)b,(double)c,(double)d));
140 :
141 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Contour points: "));
142 2 : PetscCall(RGComputeContour(rg,NPOINTS,cr,ci));
143 22 : for (i=0;i<NPOINTS;i++) {
144 : #if defined(PETSC_USE_COMPLEX)
145 : re = PetscRealPart(cr[i]);
146 : im = PetscImaginaryPart(cr[i]);
147 : #else
148 20 : re = cr[i];
149 20 : im = ci[i];
150 : #endif
151 20 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"(%.3g,%.3g) ",(double)re,(double)im));
152 : }
153 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
154 :
155 : /* check vertices */
156 2 : PetscCall(RGPolygonGetVertices(rg,&nv,&pr,&pi));
157 2 : PetscCheck(nv==NVERTEX,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Wrong number of vertices: %" PetscInt_FMT,nv);
158 16 : for (i=0;i<nv;i++) {
159 : #if !defined(PETSC_USE_COMPLEX)
160 14 : if (pr[i]!=vr[i] || pi[i]!=vi[i])
161 : #else
162 : if (pr[i]!=vr[i])
163 : #endif
164 14 : SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Vertex number %" PetscInt_FMT " does not match",i);
165 : }
166 :
167 2 : PetscCall(PetscFree(pr));
168 : #if !defined(PETSC_USE_COMPLEX)
169 2 : PetscCall(PetscFree(pi));
170 : #endif
171 2 : PetscCall(RGDestroy(&rg));
172 2 : PetscCall(SlepcFinalize());
173 : return 0;
174 : }
175 :
176 : /*TEST
177 :
178 : test:
179 : suffix: 1
180 : requires: !complex
181 :
182 : test:
183 : suffix: 1_complex
184 : requires: complex
185 :
186 : test:
187 : suffix: 2
188 : requires: !complex
189 : args: -rg_ellipse_view draw:tikz:ellipse.tikz -rg_interval_view draw:tikz:interval.tikz -rg_polygon_view draw:tikz:polygon.tikz
190 : filter: cat - ellipse.tikz interval.tikz polygon.tikz
191 : requires: !single
192 :
193 : test:
194 : suffix: 2_complex
195 : requires: complex !single
196 : args: -rg_ellipse_view draw:tikz:ellipse.tikz -rg_interval_view draw:tikz:interval.tikz -rg_polygon_view draw:tikz:polygon.tikz
197 : filter: cat - ellipse.tikz interval.tikz polygon.tikz
198 :
199 : TEST*/
|