Actual source code: rgbasic.c

slepc-main 2024-12-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: */
 10: /*
 11:    Basic RG routines
 12: */

 14: #include <slepc/private/rgimpl.h>

 16: PetscFunctionList RGList = NULL;
 17: PetscBool         RGRegisterAllCalled = PETSC_FALSE;
 18: PetscClassId      RG_CLASSID = 0;
 19: static PetscBool  RGPackageInitialized = PETSC_FALSE;

 21: /*@C
 22:    RGFinalizePackage - This function destroys everything in the Slepc interface
 23:    to the RG package. It is called from SlepcFinalize().

 25:    Level: developer

 27: .seealso: SlepcFinalize()
 28: @*/
 29: PetscErrorCode RGFinalizePackage(void)
 30: {
 31:   PetscFunctionBegin;
 32:   PetscCall(PetscFunctionListDestroy(&RGList));
 33:   RGPackageInitialized = PETSC_FALSE;
 34:   RGRegisterAllCalled  = PETSC_FALSE;
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }

 38: /*@C
 39:   RGInitializePackage - This function initializes everything in the RG package.
 40:   It is called from PetscDLLibraryRegister() when using dynamic libraries, and
 41:   on the first call to RGCreate() when using static libraries.

 43:   Level: developer

 45: .seealso: SlepcInitialize()
 46: @*/
 47: PetscErrorCode RGInitializePackage(void)
 48: {
 49:   char           logList[256];
 50:   PetscBool      opt,pkg;
 51:   PetscClassId   classids[1];

 53:   PetscFunctionBegin;
 54:   if (RGPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
 55:   RGPackageInitialized = PETSC_TRUE;
 56:   /* Register Classes */
 57:   PetscCall(PetscClassIdRegister("Region",&RG_CLASSID));
 58:   /* Register Constructors */
 59:   PetscCall(RGRegisterAll());
 60:   /* Process Info */
 61:   classids[0] = RG_CLASSID;
 62:   PetscCall(PetscInfoProcessClass("rg",1,&classids[0]));
 63:   /* Process summary exclusions */
 64:   PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt));
 65:   if (opt) {
 66:     PetscCall(PetscStrInList("rg",logList,',',&pkg));
 67:     if (pkg) PetscCall(PetscLogEventDeactivateClass(RG_CLASSID));
 68:   }
 69:   /* Register package finalizer */
 70:   PetscCall(PetscRegisterFinalize(RGFinalizePackage));
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }

 74: /*@
 75:    RGCreate - Creates an RG context.

 77:    Collective

 79:    Input Parameter:
 80: .  comm - MPI communicator

 82:    Output Parameter:
 83: .  newrg - location to put the RG context

 85:    Level: beginner

 87: .seealso: RGDestroy(), RG
 88: @*/
 89: PetscErrorCode RGCreate(MPI_Comm comm,RG *newrg)
 90: {
 91:   RG             rg;

 93:   PetscFunctionBegin;
 94:   PetscAssertPointer(newrg,2);
 95:   PetscCall(RGInitializePackage());
 96:   PetscCall(SlepcHeaderCreate(rg,RG_CLASSID,"RG","Region","RG",comm,RGDestroy,RGView));
 97:   rg->complement = PETSC_FALSE;
 98:   rg->sfactor    = 1.0;
 99:   rg->osfactor   = 0.0;
100:   rg->data       = NULL;

102:   *newrg = rg;
103:   PetscFunctionReturn(PETSC_SUCCESS);
104: }

106: /*@
107:    RGSetOptionsPrefix - Sets the prefix used for searching for all
108:    RG options in the database.

110:    Logically Collective

112:    Input Parameters:
113: +  rg     - the region context
114: -  prefix - the prefix string to prepend to all RG option requests

116:    Notes:
117:    A hyphen (-) must NOT be given at the beginning of the prefix name.
118:    The first character of all runtime options is AUTOMATICALLY the
119:    hyphen.

121:    Level: advanced

123: .seealso: RGAppendOptionsPrefix()
124: @*/
125: PetscErrorCode RGSetOptionsPrefix(RG rg,const char *prefix)
126: {
127:   PetscFunctionBegin;
129:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)rg,prefix));
130:   PetscFunctionReturn(PETSC_SUCCESS);
131: }

133: /*@
134:    RGAppendOptionsPrefix - Appends to the prefix used for searching for all
135:    RG options in the database.

137:    Logically Collective

139:    Input Parameters:
140: +  rg     - the region context
141: -  prefix - the prefix string to prepend to all RG option requests

143:    Notes:
144:    A hyphen (-) must NOT be given at the beginning of the prefix name.
145:    The first character of all runtime options is AUTOMATICALLY the hyphen.

147:    Level: advanced

149: .seealso: RGSetOptionsPrefix()
150: @*/
151: PetscErrorCode RGAppendOptionsPrefix(RG rg,const char *prefix)
152: {
153:   PetscFunctionBegin;
155:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)rg,prefix));
156:   PetscFunctionReturn(PETSC_SUCCESS);
157: }

159: /*@
160:    RGGetOptionsPrefix - Gets the prefix used for searching for all
161:    RG options in the database.

163:    Not Collective

165:    Input Parameters:
166: .  rg - the region context

168:    Output Parameters:
169: .  prefix - pointer to the prefix string used is returned

171:    Note:
172:    On the Fortran side, the user should pass in a string 'prefix' of
173:    sufficient length to hold the prefix.

175:    Level: advanced

177: .seealso: RGSetOptionsPrefix(), RGAppendOptionsPrefix()
178: @*/
179: PetscErrorCode RGGetOptionsPrefix(RG rg,const char *prefix[])
180: {
181:   PetscFunctionBegin;
183:   PetscAssertPointer(prefix,2);
184:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)rg,prefix));
185:   PetscFunctionReturn(PETSC_SUCCESS);
186: }

188: /*@
189:    RGSetType - Selects the type for the RG object.

191:    Logically Collective

193:    Input Parameters:
194: +  rg   - the region context
195: -  type - a known type

197:    Level: intermediate

199: .seealso: RGGetType()
200: @*/
201: PetscErrorCode RGSetType(RG rg,RGType type)
202: {
203:   PetscErrorCode (*r)(RG);
204:   PetscBool      match;

206:   PetscFunctionBegin;
208:   PetscAssertPointer(type,2);

210:   PetscCall(PetscObjectTypeCompare((PetscObject)rg,type,&match));
211:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

213:   PetscCall(PetscFunctionListFind(RGList,type,&r));
214:   PetscCheck(r,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested RG type %s",type);

216:   PetscTryTypeMethod(rg,destroy);
217:   PetscCall(PetscMemzero(rg->ops,sizeof(struct _RGOps)));

219:   PetscCall(PetscObjectChangeTypeName((PetscObject)rg,type));
220:   PetscCall((*r)(rg));
221:   PetscFunctionReturn(PETSC_SUCCESS);
222: }

224: /*@
225:    RGGetType - Gets the RG type name (as a string) from the RG context.

227:    Not Collective

229:    Input Parameter:
230: .  rg - the region context

232:    Output Parameter:
233: .  type - name of the region

235:    Level: intermediate

237: .seealso: RGSetType()
238: @*/
239: PetscErrorCode RGGetType(RG rg,RGType *type)
240: {
241:   PetscFunctionBegin;
243:   PetscAssertPointer(type,2);
244:   *type = ((PetscObject)rg)->type_name;
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: /*@
249:    RGSetFromOptions - Sets RG options from the options database.

251:    Collective

253:    Input Parameters:
254: .  rg - the region context

256:    Notes:
257:    To see all options, run your program with the -help option.

259:    Level: beginner

261: .seealso: RGSetOptionsPrefix()
262: @*/
263: PetscErrorCode RGSetFromOptions(RG rg)
264: {
265:   char           type[256];
266:   PetscBool      flg;
267:   PetscReal      sfactor;

269:   PetscFunctionBegin;
271:   PetscCall(RGRegisterAll());
272:   PetscObjectOptionsBegin((PetscObject)rg);
273:     PetscCall(PetscOptionsFList("-rg_type","Region type","RGSetType",RGList,(char*)(((PetscObject)rg)->type_name?((PetscObject)rg)->type_name:RGINTERVAL),type,sizeof(type),&flg));
274:     if (flg) PetscCall(RGSetType(rg,type));
275:     else if (!((PetscObject)rg)->type_name) PetscCall(RGSetType(rg,RGINTERVAL));

277:     PetscCall(PetscOptionsBool("-rg_complement","Whether region is complemented or not","RGSetComplement",rg->complement,&rg->complement,NULL));

279:     PetscCall(PetscOptionsReal("-rg_scale","Scaling factor","RGSetScale",1.0,&sfactor,&flg));
280:     if (flg) PetscCall(RGSetScale(rg,sfactor));

282:     PetscTryTypeMethod(rg,setfromoptions,PetscOptionsObject);
283:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)rg,PetscOptionsObject));
284:   PetscOptionsEnd();
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: /*@
289:    RGView - Prints the RG data structure.

291:    Collective

293:    Input Parameters:
294: +  rg - the region context
295: -  viewer - optional visualization context

297:    Note:
298:    The available visualization contexts include
299: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
300: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
301:          output where only the first processor opens
302:          the file.  All other processors send their
303:          data to the first processor to print.

305:    The user can open an alternative visualization context with
306:    PetscViewerASCIIOpen() - output to a specified file.

308:    Level: beginner

310: .seealso: RGCreate()
311: @*/
312: PetscErrorCode RGView(RG rg,PetscViewer viewer)
313: {
314:   PetscBool      isdraw,isascii;

316:   PetscFunctionBegin;
318:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)rg),&viewer));
320:   PetscCheckSameComm(rg,1,viewer,2);
321:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
322:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
323:   if (isascii) {
324:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)rg,viewer));
325:     PetscCall(PetscViewerASCIIPushTab(viewer));
326:     PetscTryTypeMethod(rg,view,viewer);
327:     PetscCall(PetscViewerASCIIPopTab(viewer));
328:     if (rg->complement) PetscCall(PetscViewerASCIIPrintf(viewer,"  selected region is the complement of the specified one\n"));
329:     if (rg->sfactor!=1.0) PetscCall(PetscViewerASCIIPrintf(viewer,"  scaling factor = %g\n",(double)rg->sfactor));
330:   } else if (isdraw) PetscTryTypeMethod(rg,view,viewer);
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: /*@
335:    RGViewFromOptions - View from options

337:    Collective

339:    Input Parameters:
340: +  rg   - the region context
341: .  obj  - optional object
342: -  name - command line option

344:    Level: intermediate

346: .seealso: RGView(), RGCreate()
347: @*/
348: PetscErrorCode RGViewFromOptions(RG rg,PetscObject obj,const char name[])
349: {
350:   PetscFunctionBegin;
352:   PetscCall(PetscObjectViewFromOptions((PetscObject)rg,obj,name));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: /*@
357:    RGIsTrivial - Whether it is the trivial region (whole complex plane).

359:    Not Collective

361:    Input Parameter:
362: .  rg - the region context

364:    Output Parameter:
365: .  trivial - true if the region is equal to the whole complex plane, e.g.,
366:              an interval region with all four endpoints unbounded or an
367:              ellipse with infinite radius.

369:    Level: beginner

371: .seealso: RGCheckInside()
372: @*/
373: PetscErrorCode RGIsTrivial(RG rg,PetscBool *trivial)
374: {
375:   PetscFunctionBegin;
378:   PetscAssertPointer(trivial,2);
379:   *trivial = PETSC_FALSE;
380:   PetscTryTypeMethod(rg,istrivial,trivial);
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

384: /*@
385:    RGCheckInside - Determines if a set of given points are inside the region or not.

387:    Not Collective

389:    Input Parameters:
390: +  rg - the region context
391: .  n  - number of points to check
392: .  ar - array of real parts
393: -  ai - array of imaginary parts

395:    Output Parameter:
396: .  inside - array of results (1=inside, 0=on the contour, -1=outside)

398:    Note:
399:    The point a is expressed as a couple of PetscScalar variables ar,ai.
400:    If built with complex scalars, the point is supposed to be stored in ar,
401:    otherwise ar,ai contain the real and imaginary parts, respectively.

403:    If a scaling factor was set, the points are scaled before checking.

405:    Level: intermediate

407: .seealso: RGSetScale(), RGSetComplement()
408: @*/
409: PetscErrorCode RGCheckInside(RG rg,PetscInt n,PetscScalar *ar,PetscScalar *ai,PetscInt *inside)
410: {
411:   PetscReal      px,py;
412:   PetscInt       i;

414:   PetscFunctionBegin;
417:   PetscAssertPointer(ar,3);
418: #if !defined(PETSC_USE_COMPLEX)
419:   PetscAssertPointer(ai,4);
420: #endif
421:   PetscAssertPointer(inside,5);

423:   for (i=0;i<n;i++) {
424: #if defined(PETSC_USE_COMPLEX)
425:     px = PetscRealPart(ar[i]);
426:     py = PetscImaginaryPart(ar[i]);
427: #else
428:     px = ar[i];
429:     py = ai[i];
430: #endif
431:     if (PetscUnlikely(rg->sfactor != 1.0)) {
432:       px /= rg->sfactor;
433:       py /= rg->sfactor;
434:     }
435:     PetscUseTypeMethod(rg,checkinside,px,py,inside+i);
436:     if (PetscUnlikely(rg->complement)) inside[i] = -inside[i];
437:   }
438:   PetscFunctionReturn(PETSC_SUCCESS);
439: }

441: /*@
442:    RGIsAxisymmetric - Determines if the region is symmetric with respect
443:    to the real or imaginary axis.

445:    Not Collective

447:    Input Parameters:
448: +  rg       - the region context
449: -  vertical - true if symmetry must be checked against the vertical axis

451:    Output Parameter:
452: .  symm - true if the region is axisymmetric

454:    Note:
455:    If the vertical argument is true, symmetry is checked with respect to
456:    the vertical axis, otherwise with respect to the horizontal axis.

458:    Level: intermediate

460: .seealso: RGCanUseConjugates()
461: @*/
462: PetscErrorCode RGIsAxisymmetric(RG rg,PetscBool vertical,PetscBool *symm)
463: {
464:   PetscFunctionBegin;
467:   PetscAssertPointer(symm,3);
468:   *symm = PETSC_FALSE;
469:   PetscTryTypeMethod(rg,isaxisymmetric,vertical,symm);
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: }

473: /*@
474:    RGCanUseConjugates - Used in contour integral methods to determine whether
475:    half of integration points can be avoided (use their conjugates).

477:    Not Collective

479:    Input Parameters:
480: +  rg       - the region context
481: -  realmats - true if the problem matrices are real

483:    Output Parameter:
484: .  useconj  - whether it is possible to use conjugates

486:    Notes:
487:    If some integration points are the conjugates of other points, then the
488:    associated computational cost can be saved. This depends on the problem
489:    matrices being real and also the region being symmetric with respect to
490:    the horizontal axis. The result is false if using real arithmetic or
491:    in the case of a flat region (height equal to zero).

493:    Level: developer

495: .seealso: RGIsAxisymmetric()
496: @*/
497: PetscErrorCode RGCanUseConjugates(RG rg,PetscBool realmats,PetscBool *useconj)
498: {
499: #if defined(PETSC_USE_COMPLEX)
500:   PetscReal      c,d;
501:   PetscBool      isaxisymm;
502: #endif

504:   PetscFunctionBegin;
507:   PetscAssertPointer(useconj,3);
508:   *useconj = PETSC_FALSE;
509: #if defined(PETSC_USE_COMPLEX)
510:   if (realmats) {
511:     PetscCall(RGIsAxisymmetric(rg,PETSC_FALSE,&isaxisymm));
512:     if (isaxisymm) {
513:       PetscCall(RGComputeBoundingBox(rg,NULL,NULL,&c,&d));
514:       if (c!=d) *useconj = PETSC_TRUE;
515:     }
516:   }
517: #endif
518:   PetscFunctionReturn(PETSC_SUCCESS);
519: }

521: /*@
522:    RGComputeContour - Computes the coordinates of several points lying on the
523:    contour of the region.

525:    Not Collective

527:    Input Parameters:
528: +  rg - the region context
529: -  n  - number of points to compute

531:    Output Parameters:
532: +  cr - location to store real parts
533: -  ci - location to store imaginary parts

535:    Notes:
536:    In real scalars, either cr or ci can be NULL (but not both). In complex
537:    scalars, the coordinates are stored in cr, which cannot be NULL (ci is
538:    not referenced).

540:    Level: intermediate

542: .seealso: RGComputeBoundingBox()
543: @*/
544: PetscErrorCode RGComputeContour(RG rg,PetscInt n,PetscScalar cr[],PetscScalar ci[])
545: {
546:   PetscInt       i;

548:   PetscFunctionBegin;
551: #if defined(PETSC_USE_COMPLEX)
552:   PetscAssertPointer(cr,3);
553: #else
554:   PetscCheck(cr || ci,PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"cr and ci cannot be NULL simultaneously");
555: #endif
556:   PetscCheck(!rg->complement,PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Cannot compute contour of region with complement flag set");
557:   PetscUseTypeMethod(rg,computecontour,n,cr,ci);
558:   for (i=0;i<n;i++) {
559:     if (cr) cr[i] *= rg->sfactor;
560:     if (ci) ci[i] *= rg->sfactor;
561:   }
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: /*@
566:    RGComputeBoundingBox - Determines the endpoints of a rectangle in the complex plane that
567:    contains the region.

569:    Not Collective

571:    Input Parameter:
572: .  rg - the region context

574:    Output Parameters:
575: +  a - left endpoint of the bounding box in the real axis
576: .  b - right endpoint of the bounding box in the real axis
577: .  c - bottom endpoint of the bounding box in the imaginary axis
578: -  d - top endpoint of the bounding box in the imaginary axis

580:    Notes:
581:    The bounding box is defined as [a,b]x[c,d]. In regions that are not bounded (e.g. an
582:    open interval) or with the complement flag set, it makes no sense to compute a bounding
583:    box, so the return values are infinite.

585:    Level: intermediate

587: .seealso: RGComputeContour()
588: @*/
589: PetscErrorCode RGComputeBoundingBox(RG rg,PetscReal *a,PetscReal *b,PetscReal *c,PetscReal *d)
590: {
591:   PetscFunctionBegin;

595:   if (rg->complement) {  /* cannot compute bounding box */
596:     if (a) *a = -PETSC_MAX_REAL;
597:     if (b) *b =  PETSC_MAX_REAL;
598:     if (c) *c = -PETSC_MAX_REAL;
599:     if (d) *d =  PETSC_MAX_REAL;
600:   } else {
601:     PetscUseTypeMethod(rg,computebbox,a,b,c,d);
602:     if (a && *a!=-PETSC_MAX_REAL) *a *= rg->sfactor;
603:     if (b && *b!= PETSC_MAX_REAL) *b *= rg->sfactor;
604:     if (c && *c!=-PETSC_MAX_REAL) *c *= rg->sfactor;
605:     if (d && *d!= PETSC_MAX_REAL) *d *= rg->sfactor;
606:   }
607:   PetscFunctionReturn(PETSC_SUCCESS);
608: }

610: /*@
611:    RGComputeQuadrature - Computes the values of the parameters used in a
612:    quadrature rule for a contour integral around the boundary of the region.

614:    Not Collective

616:    Input Parameters:
617: +  rg   - the region context
618: .  quad - the type of quadrature
619: -  n    - number of quadrature points to compute

621:    Output Parameters:
622: +  z  - quadrature points
623: .  zn - normalized quadrature points
624: -  w  - quadrature weights

626:    Notes:
627:    In complex scalars, the values returned in z are often the same as those
628:    computed by RGComputeContour(), but this is not the case in real scalars
629:    where all output arguments are real.

631:    The computed values change for different quadrature rules (trapezoidal
632:    or Chebyshev).

634:    Level: intermediate

636: .seealso: RGComputeContour()
637: @*/
638: PetscErrorCode RGComputeQuadrature(RG rg,RGQuadRule quad,PetscInt n,PetscScalar z[],PetscScalar zn[],PetscScalar w[])
639: {
640:   PetscFunctionBegin;
643:   PetscAssertPointer(z,4);
644:   PetscAssertPointer(zn,5);
645:   PetscAssertPointer(w,6);

647:   PetscCall(RGComputeContour(rg,n,z,NULL));
648:   PetscUseTypeMethod(rg,computequadrature,quad,n,z,zn,w);
649:   PetscFunctionReturn(PETSC_SUCCESS);
650: }

652: /*@
653:    RGSetComplement - Sets a flag to indicate that the region is the complement
654:    of the specified one.

656:    Logically Collective

658:    Input Parameters:
659: +  rg  - the region context
660: -  flg - the boolean flag

662:    Options Database Key:
663: .  -rg_complement <bool> - Activate/deactivate the complementation of the region

665:    Level: intermediate

667: .seealso: RGGetComplement()
668: @*/
669: PetscErrorCode RGSetComplement(RG rg,PetscBool flg)
670: {
671:   PetscFunctionBegin;
674:   rg->complement = flg;
675:   PetscFunctionReturn(PETSC_SUCCESS);
676: }

678: /*@
679:    RGGetComplement - Gets a flag that indicates whether the region
680:    is complemented or not.

682:    Not Collective

684:    Input Parameter:
685: .  rg - the region context

687:    Output Parameter:
688: .  flg - the flag

690:    Level: intermediate

692: .seealso: RGSetComplement()
693: @*/
694: PetscErrorCode RGGetComplement(RG rg,PetscBool *flg)
695: {
696:   PetscFunctionBegin;
698:   PetscAssertPointer(flg,2);
699:   *flg = rg->complement;
700:   PetscFunctionReturn(PETSC_SUCCESS);
701: }

703: /*@
704:    RGSetScale - Sets the scaling factor to be used when checking that a
705:    point is inside the region and when computing the contour.

707:    Logically Collective

709:    Input Parameters:
710: +  rg      - the region context
711: -  sfactor - the scaling factor

713:    Options Database Key:
714: .  -rg_scale <real> - Sets the scaling factor

716:    Level: advanced

718: .seealso: RGGetScale(), RGCheckInside()
719: @*/
720: PetscErrorCode RGSetScale(RG rg,PetscReal sfactor)
721: {
722:   PetscFunctionBegin;
725:   if (sfactor == (PetscReal)PETSC_DEFAULT || sfactor == (PetscReal)PETSC_DECIDE) sfactor = 1.0;
726:   PetscCheck(sfactor>0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
727:   rg->sfactor = sfactor;
728:   PetscFunctionReturn(PETSC_SUCCESS);
729: }

731: /*@
732:    RGGetScale - Gets the scaling factor.

734:    Not Collective

736:    Input Parameter:
737: .  rg - the region context

739:    Output Parameter:
740: .  sfactor - the scaling factor

742:    Level: advanced

744: .seealso: RGSetScale()
745: @*/
746: PetscErrorCode RGGetScale(RG rg,PetscReal *sfactor)
747: {
748:   PetscFunctionBegin;
750:   PetscAssertPointer(sfactor,2);
751:   *sfactor = rg->sfactor;
752:   PetscFunctionReturn(PETSC_SUCCESS);
753: }

755: /*@
756:    RGPushScale - Sets an additional scaling factor, that will multiply the
757:    user-defined scaling factor.

759:    Logically Collective

761:    Input Parameters:
762: +  rg      - the region context
763: -  sfactor - the scaling factor

765:    Notes:
766:    The current implementation does not allow pushing several scaling factors.

768:    This is intended for internal use, for instance in polynomial eigensolvers
769:    that use parameter scaling.

771:    Level: developer

773: .seealso: RGPopScale(), RGSetScale()
774: @*/
775: PetscErrorCode RGPushScale(RG rg,PetscReal sfactor)
776: {
777:   PetscFunctionBegin;
780:   PetscCheck(sfactor>0.0,PetscObjectComm((PetscObject)rg),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of scaling factor. Must be > 0");
781:   PetscCheck(!rg->osfactor,PetscObjectComm((PetscObject)rg),PETSC_ERR_SUP,"Current implementation does not allow pushing several scaling factors");
782:   rg->osfactor = rg->sfactor;
783:   rg->sfactor *= sfactor;
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

787: /*@
788:    RGPopScale - Pops the scaling factor set with RGPushScale().

790:    Logically Collective

792:    Input Parameter:
793: .  rg - the region context

795:    Level: developer

797: .seealso: RGPushScale()
798: @*/
799: PetscErrorCode RGPopScale(RG rg)
800: {
801:   PetscFunctionBegin;
803:   PetscCheck(rg->osfactor,PetscObjectComm((PetscObject)rg),PETSC_ERR_ORDER,"Must call RGPushScale first");
804:   rg->sfactor  = rg->osfactor;
805:   rg->osfactor = 0.0;
806:   PetscFunctionReturn(PETSC_SUCCESS);
807: }

809: /*@
810:    RGDestroy - Destroys RG context that was created with RGCreate().

812:    Collective

814:    Input Parameter:
815: .  rg - the region context

817:    Level: beginner

819: .seealso: RGCreate()
820: @*/
821: PetscErrorCode RGDestroy(RG *rg)
822: {
823:   PetscFunctionBegin;
824:   if (!*rg) PetscFunctionReturn(PETSC_SUCCESS);
826:   if (--((PetscObject)*rg)->refct > 0) { *rg = NULL; PetscFunctionReturn(PETSC_SUCCESS); }
827:   PetscTryTypeMethod(*rg,destroy);
828:   PetscCall(PetscHeaderDestroy(rg));
829:   PetscFunctionReturn(PETSC_SUCCESS);
830: }

832: /*@C
833:    RGRegister - Adds a region to the RG package.

835:    Not Collective

837:    Input Parameters:
838: +  name - name of a new user-defined RG
839: -  function - routine to create context

841:    Notes:
842:    RGRegister() may be called multiple times to add several user-defined regions.

844:    Level: advanced

846: .seealso: RGRegisterAll()
847: @*/
848: PetscErrorCode RGRegister(const char *name,PetscErrorCode (*function)(RG))
849: {
850:   PetscFunctionBegin;
851:   PetscCall(RGInitializePackage());
852:   PetscCall(PetscFunctionListAdd(&RGList,name,function));
853:   PetscFunctionReturn(PETSC_SUCCESS);
854: }