1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260
|
/***************************************
$Header: /cvsroot/petscgraphics/petsc.c,v 1.12 2004/07/02 20:51:08 hazelsct Exp $
This is the petsc.c main file. It has all of the PETSc-dependent functions.
***************************************/
#include <petscda.h>
#include "config.h" /* esp. for inline */
#include "illuminator.h" /* Just to make sure the interface is "right" */
extern int num_triangles;
#undef __FUNCT__
#define __FUNCT__ "DATriangulateRange"
/*++++++++++++++++++++++++++++++++++++++
Calculate vertices of isoquant triangles in a 3-D distributed array. This
takes a PETSc DA object, does some sanity checks, calculates array sizes,
gets the local vector and array, and then calls
+latex+{\tt DATriangulateLocal()}
+html+ <tt>DATriangulateLocal()</tt>
to do the rest. Note that global array access (i.e. this function) is
necessary for using default isoquant values, since we need to be able to
calculate the maximum and minimum on the global array.
int DATriangulateRange Returns 0 or an error code.
DA theda The PETSc distributed array object.
Vec globalX PETSc global vector object associated with the DA with data we'd
like to graph.
int this Index of the field we'd like to draw.
PetscScalar *minmax Position of block corners: xmin, xmax, ymin, ymax, zmin,
zmax.
int n_quants Number of isoquant surfaces to draw (isoquant values), or
+latex+{\tt PETSC\_DECIDE}
+html+ <tt>PETSC_DECIDE</tt>
to use red, yellow, green, blue at 0.2, 0.4, 0.6 and 0.8 between the vector's
global minimum and maximum values.
PetscScalar *isoquants Array of function values at which to draw isoquants,
+latex+or {\tt PETSC\_NULL} if {\tt n\_quants=PETSC\_DECIDE}.
+html+ or <tt>PETSC_NULL</tt> if <tt>n_quants=PETSC\_DECIDE</tt>.
PetscScalar *colors Array of color R,G,B,A quads for each isoquant, or
+latex+{\tt PETSC\_NULL} if {\tt n\_quants=PETSC\_DECIDE}.
+html+ <tt>PETSC_NULL</tt> if <tt>n_quants=PETSC\_DECIDE</tt>.
int xmin Smallest grid x-coordinate to render.
int xmax Largest grid x-coordinate to render, -1 goes to full x maximum, -2
in periodic systems goes to one short of x maximum.
int ymin Smallest grid y-coordinate to render.
int ymax Largest grid y-coordinate to render, -1 goes to full y maximum, -2
in periodic systems goes to one short of y maximum.
int zmin Smallest grid z-coordinate to render.
int zmax Largest grid z-coordinate to render, -1 goes to full z maximum, -2
in periodic systems goes to one short of z maximum.
++++++++++++++++++++++++++++++++++++++*/
int DATriangulateRange
(DA theda, Vec globalX, int this, PetscScalar *minmax, int n_quants,
PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin,
int ymax, int zmin,int zmax)
{
Vec localX;
PetscScalar isoquantdefaults[4],
colordefaults[16] = { 1,0,0,.5, 1,1,0,.5, 0,1,0,.5, 0,0,1,.5 };
PetscReal themax, themin;
int ierr;
/* Default isoquants */
if (n_quants == PETSC_DECIDE) {
ierr = VecStrideMin (globalX, this, PETSC_NULL, &themin); CHKERRQ (ierr);
ierr = VecStrideMax (globalX, this, PETSC_NULL, &themax); CHKERRQ (ierr);
/* Red, yellow, green, blue at 0.2, 0.4, 0.6, 0.8, all with alpha=0.5 */
n_quants = 4;
isoquantdefaults[0] = themin + 0.2*(themax-themin);
isoquantdefaults[1] = themin + 0.4*(themax-themin);
isoquantdefaults[2] = themin + 0.6*(themax-themin);
isoquantdefaults[3] = themin + 0.8*(themax-themin);
isoquants = isoquantdefaults;
colors = colordefaults;
}
/* Get the local array of points, with ghosts */
ierr = DACreateLocalVector (theda, &localX); CHKERRQ (ierr);
ierr = DAGlobalToLocalBegin (theda, globalX, INSERT_VALUES, localX); CHKERRQ(ierr);
ierr = DAGlobalToLocalEnd (theda, globalX, INSERT_VALUES, localX); CHKERRQ (ierr);
/* Use DATriangulateLocalRange() to do the work */
ierr = DATriangulateLocalRange (theda, localX, this, minmax, n_quants,
isoquants, colors, xmin,xmax, ymin,ymax,
zmin,zmax); CHKERRQ (ierr);
ierr = VecDestroy (localX); CHKERRQ (ierr);
return 0;
}
#undef __FUNCT__
#define __FUNCT__ "DATriangulateLocalRange"
/*++++++++++++++++++++++++++++++++++++++
Calculate vertices of isoquant triangles in a 3-D distributed array. This
takes a PETSc DA object, does some sanity checks, calculates array sizes, and
then gets array and passes it to Draw3DBlock for triangulation.
int DATriangulateLocalRange Returns 0 or an error code.
DA theda The PETSc distributed array object.
Vec localX PETSc local vector object associated with the DA with data we'd
like to graph.
int this Index of the field we'd like to draw.
PetscScalar *minmax Position of block corners: xmin, xmax, ymin, ymax, zmin,
zmax.
int n_quants Number of isoquant surfaces to draw (isoquant values). Note
+latex+{\tt PETSC\_DECIDE}
+html+ <tt>PETSC_DECIDE</tt>
is not a valid option here, because it's impossible to know the global
maximum and minimum and have consistent contours without user-supplied
information.
PetscScalar *isoquants Array of function values at which to draw isoquants.
PetscScalar *colors Array of color R,G,B,A quads for each isoquant.
int xmin Smallest grid x-coordinate to render.
int xmax Largest grid x-coordinate to render, -1 goes to full x maximum, -2
in periodic systems goes to one short of x maximum.
int ymin Smallest grid y-coordinate to render.
int ymax Largest grid y-coordinate to render, -1 goes to full y maximum, -2
in periodic systems goes to one short of y maximum.
int zmin Smallest grid z-coordinate to render.
int zmax Largest grid z-coordinate to render, -1 goes to full z maximum, -2
in periodic systems goes to one short of z maximum.
++++++++++++++++++++++++++++++++++++++*/
int DATriangulateLocalRange
(DA theda, Vec localX, int this, PetscScalar *minmax, int n_quants,
PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin,
int ymax, int zmin,int zmax)
{
PetscScalar *x, isoquantdefaults[4], localminmax[6],
colordefaults[16] = { 1,0,0,.5, 1,1,0,.5, 0,1,0,.5, 0,0,1,.5 };
DAStencilType stencil;
int i, ierr, fields, xw,yw,zw, xs,ys,zs, xm,ym,zm, gxs,gys,gzs, gxm,gym,gzm;
/* Default isoquant error */
if (n_quants == PETSC_DECIDE)
SETERRQ (PETSC_ERR_ARG_WRONGSTATE, "Can't get default isoquants for local array");
/* Get global and local grid boundaries */
ierr = DAGetInfo (theda, &i, &xw,&yw,&zw, PETSC_NULL,PETSC_NULL,PETSC_NULL,
&fields, PETSC_NULL, PETSC_NULL, &stencil);CHKERRQ(ierr);
if (i!=3)
SETERRQ (PETSC_ERR_ARG_WRONGSTATE, "DA must be 3-dimensional");
if (stencil!=DA_STENCIL_BOX)
SETERRQ (PETSC_ERR_ARG_WRONGSTATE, "DA must have a box stencil");
ierr = DAGetCorners (theda, &xs,&ys,&zs, &xm,&ym,&zm); CHKERRQ (ierr);
ierr = DAGetGhostCorners (theda, &gxs,&gys,&gzs, &gxm,&gym,&gzm);
CHKERRQ (ierr);
/* Get the local array of points, with ghosts */
ierr = VecGetArray (localX, &x); CHKERRQ (ierr);
/* If the array is too big, cut it down to size */
if (gxm <= xs-gxs+xm)
xm = gxm-xs+gxs-1;
if (gym <= ys-gys+ym)
ym = gym-ys+gys-1;
if (gzm <= zs-gzs+zm)
zm = gzm-zs+gzs-1;
/* Eliminate final rows/planes if *cut and periodic. */
if (xmax == -2 && xs+xm>=xw)
xm--;
if (ymax == -2 && ys+ym>=yw)
ym--;
if (zmax == -2 && zs+zm>=zw)
zm--;
/* Reframe to range */
xs = PetscMax (xs, xmin);
ys = PetscMax (ys, xmin);
zs = PetscMax (zs, xmin);
xm = (xmax > 0) ? PetscMin (xm, xmax-xs) : xm;
ym = (ymax > 0) ? PetscMin (ym, ymax-ys) : ym;
zm = (zmax > 0) ? PetscMin (zm, zmax-zs) : zm;
/* Calculate local physical size */
localminmax[0] = minmax[0] + (minmax[1]-minmax[0])*xs/xw;
localminmax[1] = minmax[2] + (minmax[1]-minmax[0])*(xs+xm)/xw;
localminmax[2] = minmax[2] + (minmax[3]-minmax[2])*ys/yw;
localminmax[3] = minmax[2] + (minmax[3]-minmax[2])*(ys+ym)/yw;
localminmax[4] = minmax[4] + (minmax[5]-minmax[4])*zs/zw;
localminmax[5] = minmax[4] + (minmax[5]-minmax[4])*(zs+zm)/zw;
/* Let 'er rip! */
ierr = Draw3DBlock (gxm,gym,gzm, xs-gxs,ys-gys,zs-gzs, xm,ym,zm, localminmax,
x+this, fields, n_quants, isoquants, colors);
CHKERRQ (ierr);
ierr = VecRestoreArray (localX, &x); CHKERRQ (ierr);
/* Prints the number of triangles on all CPUs */
#ifdef DEBUG
{
int rank;
ierr = MPI_Comm_rank (PETSC_COMM_WORLD, &rank); CHKERRQ (ierr);
ierr = PetscSynchronizedPrintf
(PETSC_COMM_WORLD, "[%d] zs=%d, zm=%d, zmin=%g, zmax=%g\n",
rank, zs, zm, localminmax[4], localminmax[5]); CHKERRQ (ierr);
ierr = PetscSynchronizedFlush (PETSC_COMM_WORLD); CHKERRQ (ierr);
ierr = PetscSynchronizedPrintf (PETSC_COMM_WORLD, "[%d] Triangles: %d\n",
rank, num_triangles); CHKERRQ (ierr);
ierr = PetscSynchronizedFlush (PETSC_COMM_WORLD); CHKERRQ (ierr);
}
#endif
return 0;
}
#undef __FUNCT__
#define __FUNCT__ "IllErrorHandler"
/*++++++++++++++++++++++++++++++++++++++
Handle errors, in this case the PETSc way.
int IllErrorHandler Returns the error code supplied.
int id Index of the error, defined in petscerror.h.
char *message Text of the error message.
++++++++++++++++++++++++++++++++++++++*/
int IllErrorHandler (int id, char *message)
{
SETERRQ (id, message);
exit (1);
}
|