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