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 the ring region.\n\n";
12 :
13 : #include <slepcrg.h>
14 :
15 : #define NPOINTS 11
16 :
17 10 : PetscErrorCode CheckPoint(RG rg,PetscReal re,PetscReal im)
18 : {
19 10 : PetscInt inside;
20 10 : PetscScalar ar,ai;
21 :
22 10 : PetscFunctionBeginUser;
23 : #if defined(PETSC_USE_COMPLEX)
24 10 : ar = PetscCMPLX(re,im);
25 : #else
26 : ar = re; ai = im;
27 : #endif
28 10 : PetscCall(RGCheckInside(rg,1,&ar,&ai,&inside));
29 16 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Point (%g,%g) is %s the region\n",(double)re,(double)im,(inside>=0)?"inside":"outside"));
30 10 : PetscFunctionReturn(PETSC_SUCCESS);
31 : }
32 :
33 5 : int main(int argc,char **argv)
34 : {
35 5 : RG rg;
36 5 : RGType rtype;
37 5 : PetscInt i;
38 5 : PetscBool triv;
39 5 : PetscReal re,im,radius,vscale,start_ang,end_ang,width,a,b,c,d;
40 5 : PetscScalar center,cr[NPOINTS],ci[NPOINTS];
41 :
42 5 : PetscFunctionBeginUser;
43 5 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
44 5 : PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
45 :
46 5 : PetscCall(RGSetType(rg,RGRING));
47 5 : PetscCall(RGIsTrivial(rg,&triv));
48 5 : PetscCheck(triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be trivial before setting parameters");
49 5 : PetscCall(RGRingSetParameters(rg,2,PETSC_DETERMINE,0.5,0.25,0.75,0.1));
50 5 : PetscCall(RGSetFromOptions(rg));
51 5 : PetscCall(RGIsTrivial(rg,&triv));
52 5 : PetscCheck(!triv,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Region should be non-trivial after setting parameters");
53 5 : PetscCall(RGView(rg,NULL));
54 5 : PetscCall(RGViewFromOptions(rg,NULL,"-rg_view"));
55 :
56 5 : PetscCall(RGGetType(rg,&rtype));
57 5 : PetscCall(RGRingGetParameters(rg,¢er,&radius,&vscale,&start_ang,&end_ang,&width));
58 5 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s region: \n center=%g, radius=%g, vscale=%g\n start angle=%g, end angle=%g, width=%g\n\n",rtype,(double)PetscRealPart(center),(double)radius,(double)vscale,(double)start_ang,(double)end_ang,(double)width));
59 :
60 5 : PetscCall(CheckPoint(rg,3.0,0.3));
61 5 : PetscCall(CheckPoint(rg,1.1747,0.28253));
62 :
63 5 : PetscCall(RGComputeBoundingBox(rg,&a,&b,&c,&d));
64 5 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"The bounding box is [%g,%g]x[%g,%g]\n",(double)a,(double)b,(double)c,(double)d));
65 :
66 5 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Contour points: "));
67 5 : PetscCall(RGComputeContour(rg,NPOINTS,cr,ci));
68 60 : for (i=0;i<NPOINTS;i++) {
69 : #if defined(PETSC_USE_COMPLEX)
70 55 : re = PetscRealPart(cr[i]);
71 55 : im = PetscImaginaryPart(cr[i]);
72 : #else
73 : re = cr[i];
74 : im = ci[i];
75 : #endif
76 55 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"(%.3g,%.3g) ",(double)re,(double)im));
77 : }
78 5 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
79 :
80 5 : PetscCall(RGDestroy(&rg));
81 5 : PetscCall(SlepcFinalize());
82 : return 0;
83 : }
84 :
85 : /*TEST
86 :
87 : test:
88 : suffix: 1
89 : args: -rg_ring_width 0.015
90 :
91 : test:
92 : suffix: 2
93 : args: -rg_ring_width 0.015 -rg_scale 1.5
94 :
95 : test:
96 : suffix: 3
97 : args: -rg_view draw:tikz:test2_3_ring.tikz
98 : filter: cat - test2_3_ring.tikz
99 : requires: !single
100 :
101 : test:
102 : suffix: 4
103 : args: -rg_ring_width 0.015 -rg_ring_center 3 -rg_ring_radius 0.3 -rg_ring_vscale 1
104 : requires: !single
105 :
106 : test:
107 : suffix: 5
108 : args: -rg_ring_width 0.1 -rg_ring_center 0.35 -rg_ring_radius 0.86 -rg_ring_vscale 1 -rg_ring_startangle 0.75 -rg_ring_endangle 0.25
109 :
110 : TEST*/
|