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 : Ring region, similar to the ellipse but with a start and end angle,
12 : together with the width
13 : */
14 :
15 : #include <slepc/private/rgimpl.h> /*I "slepcrg.h" I*/
16 : #include <petscdraw.h>
17 :
18 : typedef struct {
19 : PetscScalar center; /* center of the ellipse */
20 : PetscReal radius; /* radius of the ellipse */
21 : PetscReal vscale; /* vertical scale of the ellipse */
22 : PetscReal start_ang; /* start angle */
23 : PetscReal end_ang; /* end angle */
24 : PetscReal width; /* ring width */
25 : } RG_RING;
26 :
27 19 : static PetscErrorCode RGRingSetParameters_Ring(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale,PetscReal start_ang,PetscReal end_ang,PetscReal width)
28 : {
29 19 : RG_RING *ctx = (RG_RING*)rg->data;
30 :
31 19 : PetscFunctionBegin;
32 19 : ctx->center = center;
33 19 : if (radius == (PetscReal)PETSC_DEFAULT) {
34 5 : ctx->radius = 1.0;
35 : } else {
36 14 : PetscCheck(radius>0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The radius argument must be > 0.0");
37 14 : ctx->radius = radius;
38 : }
39 19 : PetscCheck(vscale>0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The vscale argument must be > 0.0");
40 19 : ctx->vscale = vscale;
41 19 : if (start_ang == (PetscReal)PETSC_DEFAULT) {
42 0 : ctx->start_ang = 0.0;
43 : } else {
44 19 : PetscCheck(start_ang>=0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The right-hand side angle argument must be >= 0.0");
45 19 : PetscCheck(start_ang<=1.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The right-hand side angle argument must be <= 1.0");
46 19 : ctx->start_ang = start_ang;
47 : }
48 19 : if (end_ang == (PetscReal)PETSC_DEFAULT) {
49 0 : ctx->end_ang = 1.0;
50 : } else {
51 19 : PetscCheck(end_ang>=0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The left-hand side angle argument must be >= 0.0");
52 19 : PetscCheck(end_ang<=1.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The left-hand side angle argument must be <= 1.0");
53 19 : ctx->end_ang = end_ang;
54 : }
55 : #if !defined(PETSC_USE_COMPLEX)
56 : PetscCheck(ctx->start_ang+ctx->end_ang==1.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_WRONG,"In real scalars the region must be symmetric wrt real axis");
57 : #endif
58 19 : if (width == (PetscReal)PETSC_DEFAULT) {
59 0 : ctx->width = 0.1;
60 : } else {
61 19 : PetscCheck(width>0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"The width argument must be > 0.0");
62 19 : ctx->width = width;
63 : }
64 19 : PetscFunctionReturn(PETSC_SUCCESS);
65 : }
66 :
67 : /*@
68 : RGRingSetParameters - Sets the parameters defining the ring region.
69 :
70 : Logically Collective
71 :
72 : Input Parameters:
73 : + rg - the region context
74 : . center - center of the ellipse
75 : . radius - radius of the ellipse
76 : . vscale - vertical scale of the ellipse
77 : . start_ang - the right-hand side angle
78 : . end_ang - the left-hand side angle
79 : - width - width of the ring
80 :
81 : Options Database Keys:
82 : + -rg_ring_center - Sets the center
83 : . -rg_ring_radius - Sets the radius
84 : . -rg_ring_vscale - Sets the vertical scale
85 : . -rg_ring_startangle - Sets the right-hand side angle
86 : . -rg_ring_endangle - Sets the left-hand side angle
87 : - -rg_ring_width - Sets the width of the ring
88 :
89 : Notes:
90 : The values of center, radius and vscale have the same meaning as in the
91 : ellipse region. The startangle and endangle define the span of the ring
92 : (by default it is the whole ring), while the width is the separation
93 : between the two concentric ellipses (above and below the radius by
94 : width/2).
95 :
96 : The start and end angles are expressed as a fraction of the circumference.
97 : The allowed range is [0..1], with 0 corresponding to 0 radians, 0.25 to
98 : pi/2 radians, and so on. It is allowed to have startangle>endangle, in
99 : which case the ring region crosses over the zero angle.
100 :
101 : In the case of complex scalars, a complex center can be provided in the
102 : command line with [+/-][realnumber][+/-]realnumberi with no spaces, e.g.
103 : -rg_ring_center 1.0+2.0i
104 :
105 : When PETSc is built with real scalars, the center is restricted to a real value,
106 : and the start and end angles must be such that the region is symmetric with
107 : respect to the real axis.
108 :
109 : Level: advanced
110 :
111 : .seealso: RGRingGetParameters()
112 : @*/
113 19 : PetscErrorCode RGRingSetParameters(RG rg,PetscScalar center,PetscReal radius,PetscReal vscale,PetscReal start_ang,PetscReal end_ang,PetscReal width)
114 : {
115 19 : PetscFunctionBegin;
116 19 : PetscValidHeaderSpecific(rg,RG_CLASSID,1);
117 76 : PetscValidLogicalCollectiveScalar(rg,center,2);
118 76 : PetscValidLogicalCollectiveReal(rg,radius,3);
119 76 : PetscValidLogicalCollectiveReal(rg,vscale,4);
120 76 : PetscValidLogicalCollectiveReal(rg,start_ang,5);
121 76 : PetscValidLogicalCollectiveReal(rg,end_ang,6);
122 76 : PetscValidLogicalCollectiveReal(rg,width,7);
123 19 : PetscTryMethod(rg,"RGRingSetParameters_C",(RG,PetscScalar,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal),(rg,center,radius,vscale,start_ang,end_ang,width));
124 19 : PetscFunctionReturn(PETSC_SUCCESS);
125 : }
126 :
127 24 : static PetscErrorCode RGRingGetParameters_Ring(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale,PetscReal *start_ang,PetscReal *end_ang,PetscReal *width)
128 : {
129 24 : RG_RING *ctx = (RG_RING*)rg->data;
130 :
131 24 : PetscFunctionBegin;
132 24 : if (center) *center = ctx->center;
133 24 : if (radius) *radius = ctx->radius;
134 24 : if (vscale) *vscale = ctx->vscale;
135 24 : if (start_ang) *start_ang = ctx->start_ang;
136 24 : if (end_ang) *end_ang = ctx->end_ang;
137 24 : if (width) *width = ctx->width;
138 24 : PetscFunctionReturn(PETSC_SUCCESS);
139 : }
140 :
141 : /*@
142 : RGRingGetParameters - Gets the parameters that define the ring region.
143 :
144 : Not Collective
145 :
146 : Input Parameter:
147 : . rg - the region context
148 :
149 : Output Parameters:
150 : + center - center of the region
151 : . radius - radius of the region
152 : . vscale - vertical scale of the region
153 : . start_ang - the right-hand side angle
154 : . end_ang - the left-hand side angle
155 : - width - width of the ring
156 :
157 : Level: advanced
158 :
159 : .seealso: RGRingSetParameters()
160 : @*/
161 24 : PetscErrorCode RGRingGetParameters(RG rg,PetscScalar *center,PetscReal *radius,PetscReal *vscale,PetscReal *start_ang,PetscReal *end_ang,PetscReal *width)
162 : {
163 24 : PetscFunctionBegin;
164 24 : PetscValidHeaderSpecific(rg,RG_CLASSID,1);
165 24 : PetscUseMethod(rg,"RGRingGetParameters_C",(RG,PetscScalar*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(rg,center,radius,vscale,start_ang,end_ang,width));
166 24 : PetscFunctionReturn(PETSC_SUCCESS);
167 : }
168 :
169 6 : static PetscErrorCode RGView_Ring(RG rg,PetscViewer viewer)
170 : {
171 6 : RG_RING *ctx = (RG_RING*)rg->data;
172 6 : int winw,winh;
173 6 : PetscBool isdraw,isascii;
174 6 : PetscDraw draw;
175 6 : PetscDrawAxis axis;
176 6 : PetscReal cx,cy,radius,width,ab,cd,lx,ly,w,end_ang,x1,y1,x2,y2,r,theta,scale=1.2;
177 6 : char str[50];
178 :
179 6 : PetscFunctionBegin;
180 6 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
181 6 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
182 6 : if (isascii) {
183 5 : PetscCall(SlepcSNPrintfScalar(str,sizeof(str),ctx->center,PETSC_FALSE));
184 5 : PetscCall(PetscViewerASCIIPrintf(viewer," center: %s, radius: %g, vscale: %g, start angle: %g, end angle: %g, ring width: %g\n",str,RGShowReal(ctx->radius),RGShowReal(ctx->vscale),(double)ctx->start_ang,(double)ctx->end_ang,(double)ctx->width));
185 1 : } else if (isdraw) {
186 1 : PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
187 1 : PetscCall(PetscDrawCheckResizedWindow(draw));
188 1 : PetscCall(PetscDrawGetWindowSize(draw,&winw,&winh));
189 1 : winw = PetscMax(winw,1); winh = PetscMax(winh,1);
190 1 : PetscCall(PetscDrawClear(draw));
191 1 : PetscCall(PetscDrawSetTitle(draw,"Ring region"));
192 1 : PetscCall(PetscDrawAxisCreate(draw,&axis));
193 1 : cx = PetscRealPart(ctx->center)*rg->sfactor;
194 1 : cy = PetscImaginaryPart(ctx->center)*rg->sfactor;
195 1 : radius = ctx->radius*rg->sfactor;
196 1 : width = ctx->width*rg->sfactor;
197 1 : lx = 2*(radius+width);
198 1 : ly = 2*(radius+width)*ctx->vscale;
199 1 : ab = cx;
200 1 : cd = cy;
201 1 : w = scale*PetscMax(lx/winw,ly/winh)/2;
202 1 : PetscCall(PetscDrawAxisSetLimits(axis,ab-w*winw,ab+w*winw,cd-w*winh,cd+w*winh));
203 1 : PetscCall(PetscDrawAxisDraw(axis));
204 1 : PetscCall(PetscDrawAxisDestroy(&axis));
205 : /* draw outer ellipse */
206 1 : PetscCall(PetscDrawEllipse(draw,cx,cy,2*(radius+width),2*(radius+width)*ctx->vscale,PETSC_DRAW_ORANGE));
207 : /* remove inner part */
208 1 : PetscCall(PetscDrawEllipse(draw,cx,cy,2*(radius-width),2*(radius-width)*ctx->vscale,PETSC_DRAW_WHITE));
209 1 : if (ctx->start_ang!=ctx->end_ang) {
210 : /* remove section from end_ang to start_ang */
211 1 : end_ang = (ctx->start_ang<ctx->end_ang)? ctx->end_ang-1: ctx->end_ang;
212 1 : theta = end_ang;
213 1 : r = scale*(radius+width);
214 1 : if (ctx->vscale>1) r *= ctx->vscale;
215 2 : x1 = PetscMin(PetscMax(ab+r*PetscCosReal(2.0*PETSC_PI*theta),ab-w*winw),ab+w*winw);
216 2 : y1 = PetscMin(PetscMax(cd+r*PetscSinReal(2.0*PETSC_PI*theta),cd-w*winh),cd+w*winh);
217 4 : do {
218 4 : theta = PetscMin(PetscFloorReal(8*theta+1)/8,ctx->start_ang);
219 8 : x2 = PetscMin(PetscMax(ab+r*PetscCosReal(2.0*PETSC_PI*theta),ab-w*winw),ab+w*winw);
220 8 : y2 = PetscMin(PetscMax(cd+r*PetscSinReal(2.0*PETSC_PI*theta),cd-w*winh),cd+w*winh);
221 4 : PetscCall(PetscDrawTriangle(draw,cx,cy,x1,y1,x2,y2,PETSC_DRAW_WHITE,PETSC_DRAW_WHITE,PETSC_DRAW_WHITE));
222 4 : x1 = x2; y1 = y2;
223 4 : } while (theta<ctx->start_ang);
224 : }
225 1 : PetscCall(PetscDrawFlush(draw));
226 1 : PetscCall(PetscDrawSave(draw));
227 1 : PetscCall(PetscDrawPause(draw));
228 : }
229 6 : PetscFunctionReturn(PETSC_SUCCESS);
230 : }
231 :
232 13 : static PetscErrorCode RGIsTrivial_Ring(RG rg,PetscBool *trivial)
233 : {
234 13 : RG_RING *ctx = (RG_RING*)rg->data;
235 :
236 13 : PetscFunctionBegin;
237 13 : if (rg->complement) *trivial = PetscNot(ctx->radius);
238 13 : else *trivial = PetscNot(ctx->radius<PETSC_MAX_REAL);
239 13 : PetscFunctionReturn(PETSC_SUCCESS);
240 : }
241 :
242 8 : static PetscErrorCode RGComputeContour_Ring(RG rg,PetscInt n,PetscScalar *cr,PetscScalar *ci)
243 : {
244 8 : RG_RING *ctx = (RG_RING*)rg->data;
245 8 : PetscReal theta,start_ang;
246 8 : PetscInt i,n2=n/2;
247 :
248 8 : PetscFunctionBegin;
249 8 : start_ang = (ctx->start_ang>ctx->end_ang)? ctx->start_ang-1: ctx->start_ang;
250 159 : for (i=0;i<n;i++) {
251 151 : if (i < n2) {
252 73 : theta = ((ctx->end_ang-start_ang)*i/n2 + start_ang)*2.0*PETSC_PI;
253 : #if defined(PETSC_USE_COMPLEX)
254 73 : cr[i] = ctx->center + (ctx->radius+ctx->width/2.0)*PetscCMPLX(PetscCosReal(theta),ctx->vscale*PetscSinReal(theta));
255 : #else
256 : if (cr) cr[i] = ctx->center + (ctx->radius+ctx->width/2.0)*PetscCosReal(theta);
257 : if (ci) ci[i] = (ctx->radius+ctx->width/2.0)*ctx->vscale*PetscSinReal(theta);
258 : #endif
259 : } else {
260 78 : theta = ((ctx->end_ang-start_ang)*(n-i)/n2 + start_ang)*2.0*PETSC_PI;
261 : #if defined(PETSC_USE_COMPLEX)
262 78 : cr[i] = ctx->center + (ctx->radius-ctx->width/2.0)*PetscCMPLX(PetscCosReal(theta),ctx->vscale*PetscSinReal(theta));
263 : #else
264 : if (cr) cr[i] = ctx->center + (ctx->radius-ctx->width/2.0)*PetscCosReal(theta);
265 : if (ci) ci[i] = (ctx->radius-ctx->width/2.0)*ctx->vscale*PetscSinReal(theta);
266 : #endif
267 : }
268 : }
269 8 : PetscFunctionReturn(PETSC_SUCCESS);
270 : }
271 :
272 5 : static PetscErrorCode RGComputeBoundingBox_Ring(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
273 : {
274 5 : RG_RING *ctx = (RG_RING*)rg->data;
275 :
276 5 : PetscFunctionBegin;
277 : /* current implementation does not return a tight bounding box */
278 5 : if (a) *a = PetscRealPart(ctx->center) - (ctx->radius+ctx->width/2.0);
279 5 : if (b) *b = PetscRealPart(ctx->center) + (ctx->radius+ctx->width/2.0);
280 5 : if (c) *c = PetscImaginaryPart(ctx->center) - (ctx->radius+ctx->width/2.0)*ctx->vscale;
281 5 : if (d) *d = PetscImaginaryPart(ctx->center) + (ctx->radius+ctx->width/2.0)*ctx->vscale;
282 5 : PetscFunctionReturn(PETSC_SUCCESS);
283 : }
284 :
285 3 : static PetscErrorCode RGComputeQuadrature_Ring(RG rg,RGQuadRule quad,PetscInt n,PetscScalar *z,PetscScalar *zn,PetscScalar *w)
286 : {
287 3 : RG_RING *ctx = (RG_RING*)rg->data;
288 3 : PetscReal max_w=0.0;
289 3 : PetscScalar tmp,tmp2;
290 3 : PetscInt i,j;
291 :
292 3 : PetscFunctionBegin;
293 3 : if (quad == RG_QUADRULE_CHEBYSHEV) {
294 : #if defined(PETSC_USE_COMPLEX)
295 : PetscReal theta;
296 33 : for (i=0;i<n;i++) {
297 32 : theta = PETSC_PI*(i+0.5)/n;
298 32 : zn[i] = PetscCosReal(theta);
299 32 : w[i] = PetscCosReal((n-1)*theta)/n;
300 32 : theta = (ctx->start_ang*2.0+(ctx->end_ang-ctx->start_ang)*(PetscRealPart(zn[i])+1.0))*PETSC_PI;
301 32 : z[i] = rg->sfactor*(ctx->center + ctx->radius*PetscCMPLX(PetscCosReal(theta),ctx->vscale*PetscSinReal(theta)));
302 : }
303 : #else
304 : SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Integration points on a vertical line require complex arithmetic");
305 : #endif
306 : } else { /* RG_QUADRULE_TRAPEZOIDAL */
307 66 : for (i=0;i<n;i++) {
308 64 : zn[i] = (z[i]-rg->sfactor*ctx->center)/(rg->sfactor*ctx->radius);
309 64 : tmp = 1.0; tmp2 = 1.0;
310 2112 : for (j=0;j<n;j++) {
311 2048 : tmp *= z[j];
312 2048 : if (i != j) tmp2 *= z[j]-z[i];
313 : }
314 64 : w[i] = tmp/tmp2;
315 86 : max_w = PetscMax(PetscAbsScalar(w[i]),max_w);
316 : }
317 66 : for (i=0;i<n;i++) w[i] /= (PetscScalar)max_w;
318 : }
319 3 : PetscFunctionReturn(PETSC_SUCCESS);
320 : }
321 :
322 69 : static PetscErrorCode RGCheckInside_Ring(RG rg,PetscReal px,PetscReal py,PetscInt *inside)
323 : {
324 69 : RG_RING *ctx = (RG_RING*)rg->data;
325 69 : PetscReal dx,dy,r;
326 :
327 69 : PetscFunctionBegin;
328 : /* outer ellipse */
329 : #if defined(PETSC_USE_COMPLEX)
330 69 : dx = (px-PetscRealPart(ctx->center))/(ctx->radius+ctx->width/2.0);
331 69 : dy = (py-PetscImaginaryPart(ctx->center))/(ctx->radius+ctx->width/2.0);
332 : #else
333 : dx = (px-ctx->center)/(ctx->radius+ctx->width/2.0);
334 : dy = py/(ctx->radius+ctx->width/2.0);
335 : #endif
336 69 : r = 1.0-dx*dx-(dy*dy)/(ctx->vscale*ctx->vscale);
337 69 : *inside = PetscSign(r);
338 : /* inner ellipse */
339 : #if defined(PETSC_USE_COMPLEX)
340 69 : dx = (px-PetscRealPart(ctx->center))/(ctx->radius-ctx->width/2.0);
341 69 : dy = (py-PetscImaginaryPart(ctx->center))/(ctx->radius-ctx->width/2.0);
342 : #else
343 : dx = (px-ctx->center)/(ctx->radius-ctx->width/2.0);
344 : dy = py/(ctx->radius-ctx->width/2.0);
345 : #endif
346 69 : r = -1.0+dx*dx+(dy*dy)/(ctx->vscale*ctx->vscale);
347 69 : *inside *= PetscSign(r);
348 69 : if (*inside == 1) { /* check angles */
349 : #if defined(PETSC_USE_COMPLEX)
350 45 : dx = (px-PetscRealPart(ctx->center));
351 45 : dy = (py-PetscImaginaryPart(ctx->center));
352 : #else
353 : dx = px-ctx->center;
354 : dy = py;
355 : #endif
356 45 : if (dx == 0) {
357 1 : if (dy == 0) r = -1;
358 1 : else if (dy > 0) r = 0.25;
359 0 : else r = 0.75;
360 44 : } else if (dx > 0) {
361 25 : r = PetscAtanReal((dy/ctx->vscale)/dx);
362 25 : if (dy >= 0) r /= 2*PETSC_PI;
363 12 : else r = r/(2*PETSC_PI)+1;
364 19 : } else r = PetscAtanReal((dy/ctx->vscale)/dx)/(2*PETSC_PI)+0.5;
365 45 : if (ctx->start_ang>ctx->end_ang) {
366 25 : if (r>ctx->end_ang && r<ctx->start_ang) *inside = -1;
367 : } else {
368 20 : if (r<ctx->start_ang || r>ctx->end_ang) *inside = -1;
369 : }
370 : }
371 69 : PetscFunctionReturn(PETSC_SUCCESS);
372 : }
373 :
374 4 : static PetscErrorCode RGIsAxisymmetric_Ring(RG rg,PetscBool vertical,PetscBool *symm)
375 : {
376 4 : RG_RING *ctx = (RG_RING*)rg->data;
377 :
378 4 : PetscFunctionBegin;
379 5 : if (vertical) *symm = (PetscRealPart(ctx->center) == 0.0 && PetscAbs(ctx->start_ang+ctx->end_ang-PetscRealConstant(1.0)) == 0.5)? PETSC_TRUE: PETSC_FALSE;
380 3 : else *symm = (PetscImaginaryPart(ctx->center) == 0.0 && ctx->start_ang+ctx->end_ang == 1.0)? PETSC_TRUE: PETSC_FALSE;
381 4 : PetscFunctionReturn(PETSC_SUCCESS);
382 : }
383 :
384 15 : static PetscErrorCode RGSetFromOptions_Ring(RG rg,PetscOptionItems *PetscOptionsObject)
385 : {
386 15 : PetscScalar s;
387 15 : PetscReal r1,r2,r3,r4,r5;
388 15 : PetscBool flg1,flg2,flg3,flg4,flg5,flg6;
389 :
390 15 : PetscFunctionBegin;
391 15 : PetscOptionsHeadBegin(PetscOptionsObject,"RG Ring Options");
392 :
393 15 : PetscCall(RGRingGetParameters(rg,&s,&r1,&r2,&r3,&r4,&r5));
394 15 : PetscCall(PetscOptionsScalar("-rg_ring_center","Center of ellipse","RGRingSetParameters",s,&s,&flg1));
395 15 : PetscCall(PetscOptionsReal("-rg_ring_radius","Radius of ellipse","RGRingSetParameters",r1,&r1,&flg2));
396 15 : PetscCall(PetscOptionsReal("-rg_ring_vscale","Vertical scale of ellipse","RGRingSetParameters",r2,&r2,&flg3));
397 15 : PetscCall(PetscOptionsReal("-rg_ring_startangle","Right-hand side angle","RGRingSetParameters",r3,&r3,&flg4));
398 15 : PetscCall(PetscOptionsReal("-rg_ring_endangle","Left-hand side angle","RGRingSetParameters",r4,&r4,&flg5));
399 15 : PetscCall(PetscOptionsReal("-rg_ring_width","Width of ring","RGRingSetParameters",r5,&r5,&flg6));
400 15 : if (flg1 || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(RGRingSetParameters(rg,s,r1,r2,r3,r4,r5));
401 :
402 15 : PetscOptionsHeadEnd();
403 15 : PetscFunctionReturn(PETSC_SUCCESS);
404 : }
405 :
406 12 : static PetscErrorCode RGDestroy_Ring(RG rg)
407 : {
408 12 : PetscFunctionBegin;
409 12 : PetscCall(PetscFree(rg->data));
410 12 : PetscCall(PetscObjectComposeFunction((PetscObject)rg,"RGRingSetParameters_C",NULL));
411 12 : PetscCall(PetscObjectComposeFunction((PetscObject)rg,"RGRingGetParameters_C",NULL));
412 12 : PetscFunctionReturn(PETSC_SUCCESS);
413 : }
414 :
415 12 : SLEPC_EXTERN PetscErrorCode RGCreate_Ring(RG rg)
416 : {
417 12 : RG_RING *ring;
418 :
419 12 : PetscFunctionBegin;
420 12 : PetscCall(PetscNew(&ring));
421 12 : ring->center = 0.0;
422 12 : ring->radius = PETSC_MAX_REAL;
423 12 : ring->vscale = 1.0;
424 12 : ring->start_ang = 0.0;
425 12 : ring->end_ang = 1.0;
426 12 : ring->width = 0.1;
427 12 : rg->data = (void*)ring;
428 :
429 12 : rg->ops->istrivial = RGIsTrivial_Ring;
430 12 : rg->ops->computecontour = RGComputeContour_Ring;
431 12 : rg->ops->computebbox = RGComputeBoundingBox_Ring;
432 12 : rg->ops->computequadrature = RGComputeQuadrature_Ring;
433 12 : rg->ops->checkinside = RGCheckInside_Ring;
434 12 : rg->ops->isaxisymmetric = RGIsAxisymmetric_Ring;
435 12 : rg->ops->setfromoptions = RGSetFromOptions_Ring;
436 12 : rg->ops->view = RGView_Ring;
437 12 : rg->ops->destroy = RGDestroy_Ring;
438 12 : PetscCall(PetscObjectComposeFunction((PetscObject)rg,"RGRingSetParameters_C",RGRingSetParameters_Ring));
439 12 : PetscCall(PetscObjectComposeFunction((PetscObject)rg,"RGRingGetParameters_C",RGRingGetParameters_Ring));
440 12 : PetscFunctionReturn(PETSC_SUCCESS);
441 : }
|