Actual source code: test2.c

slepc-3.17.0 2022-03-31
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: */

 11: static char help[] = "Test the ring region.\n\n";

 13: #include <slepcrg.h>

 15: #define NPOINTS 11

 17: PetscErrorCode CheckPoint(RG rg,PetscReal re,PetscReal im)
 18: {
 19:   PetscInt       inside;
 20:   PetscScalar    ar,ai;

 23: #if defined(PETSC_USE_COMPLEX)
 24:   ar = PetscCMPLX(re,im);
 25: #else
 26:   ar = re; ai = im;
 27: #endif
 28:   RGCheckInside(rg,1,&ar,&ai,&inside);
 29:   PetscPrintf(PETSC_COMM_WORLD,"Point (%g,%g) is %s the region\n",(double)re,(double)im,(inside>=0)?"inside":"outside");
 30:   PetscFunctionReturn(0);
 31: }

 33: int main(int argc,char **argv)
 34: {
 35:   RG             rg;
 36:   RGType         rtype;
 37:   PetscInt       i;
 38:   PetscBool      triv;
 39:   PetscReal      re,im,radius,vscale,start_ang,end_ang,width,a,b,c,d;
 40:   PetscScalar    center,cr[NPOINTS],ci[NPOINTS];

 42:   SlepcInitialize(&argc,&argv,(char*)0,help);
 43:   RGCreate(PETSC_COMM_WORLD,&rg);

 45:   RGSetType(rg,RGRING);
 46:   RGIsTrivial(rg,&triv);
 48:   RGRingSetParameters(rg,2,PETSC_DEFAULT,0.5,0.25,0.75,0.1);
 49:   RGSetFromOptions(rg);
 50:   RGIsTrivial(rg,&triv);
 52:   RGView(rg,NULL);
 53:   RGViewFromOptions(rg,NULL,"-rg_view");

 55:   RGGetType(rg,&rtype);
 56:   RGRingGetParameters(rg,&center,&radius,&vscale,&start_ang,&end_ang,&width);
 57:   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:   CheckPoint(rg,3.0,0.3);
 60:   CheckPoint(rg,1.1747,0.28253);

 62:   RGComputeBoundingBox(rg,&a,&b,&c,&d);
 63:   PetscPrintf(PETSC_COMM_WORLD,"The bounding box is [%g,%g]x[%g,%g]\n",(double)a,(double)b,(double)c,(double)d);

 65:   PetscPrintf(PETSC_COMM_WORLD,"Contour points: ");
 66:   RGComputeContour(rg,NPOINTS,cr,ci);
 67:   for (i=0;i<NPOINTS;i++) {
 68: #if defined(PETSC_USE_COMPLEX)
 69:     re = PetscRealPart(cr[i]);
 70:     im = PetscImaginaryPart(cr[i]);
 71: #else
 72:     re = cr[i];
 73:     im = ci[i];
 74: #endif
 75:     PetscPrintf(PETSC_COMM_WORLD,"(%.3g,%.3g) ",(double)re,(double)im);
 76:   }
 77:   PetscPrintf(PETSC_COMM_WORLD,"\n");

 79:   RGDestroy(&rg);
 80:   SlepcFinalize();
 81:   return 0;
 82: }

 84: /*TEST

 86:    test:
 87:       suffix: 1
 88:       args: -rg_ring_width 0.015

 90:    test:
 91:       suffix: 2
 92:       args: -rg_ring_width 0.015 -rg_scale 1.5

 94:    test:
 95:       suffix: 3
 96:       args: -rg_view draw:tikz:test2_3_ring.tikz
 97:       filter: cat - test2_3_ring.tikz
 98:       requires: !single

100:    test:
101:       suffix: 4
102:       args: -rg_ring_width 0.015 -rg_ring_center 3 -rg_ring_radius 0.3 -rg_ring_vscale 1
103:       requires: !single

105:    test:
106:       suffix: 5
107:       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: TEST*/