1    | /***************************************
2    |   $Header$
3    | 
4    |   This is the interface for the Illuminator library.
5    | ***************************************/
6    | 
7    | #ifndef ILLUMINATOR_H
8    | #define ILLUMINATOR_H    /*+ To stop multiple inclusions. +*/
9    | #include <petscda.h>
10   | #include <glib.h>
11   | #include <gtk/gtk.h>
12   | 
13   | /*+ The
14   |   +latex+{\tt ISurface}
15   |   +html+ <tt>ISurface</tt>
16   |   type is the container (or object class) for triangle data which illuminator
17   |   will render. +*/
18   | typedef void * ISurface;
19   | 
20   | int ISurfCreate (ISurface *newsurf);
21   | int ISurfDestroy (ISurface Surf);
22   | int ISurfClear (ISurface Surf);
23   | 
24   | /*+ The
25   |   +latex+{\tt IDisplay}
26   |   +html+ <tt>IDisplay</tt>
27   |   type is the container for display data, such as the geomview output pipe, RGB
28   |   buffer, multi-layer z-buffer, etc. +*/
29   | typedef void * IDisplay;
30   | 
31   | int IDispCreate (IDisplay *newdisp, int width,int height,int rowskip,int bpp,
32   | 		 int zbuf);
33   | int IDispResize (IDisplay Disp, int width,int height,int rowskip,int bpp,
34   | 		 int zbuf);
35   | int IDispDestroy (IDisplay Disp);
36   | int IDispFill (IDisplay Disp, guchar *color);
37   | int IDispDrawGdk (IDisplay Disp, GtkWidget *dataview, GdkRgbDither dith);
38   | int IDispWritePPM (IDisplay Disp, char *filename);
39   | 
40   | /*+ A value of
41   |   +latex+{\tt field\_plot\_type}
42   |   +html+ <tt>field_plot_type</tt>
43   |   is attached to each field in a simulation in order to visualize them
44   |   properly.  Types are as follows:
45   |   +*/
46   | typedef enum {
47   |   /*+Scalar field.+*/
48   |   FIELD_SCALAR           = 0x00,
49   |   /*+Ternary composition field with two components (third component is inferred
50   |     from first two).+*/
51   |   FIELD_TERNARY          = 0x10,
52   |   /*+Ternary composition with pseudo-components mapping onto a rectangle
53   |     instead of a triangle.+*/
54   |   FIELD_TERNARY_SQUARE   = 0x18,
55   |   /*+Vector field.+*/
56   |   FIELD_VECTOR           = 0x20,
57   |   /*+Full ds*ds tensor field, e.g. transformation.+*/
58   |   FIELD_TENSOR_FULL      = 0x30,
59   |   /*+Symmetric tensor field (using lines in principal stress directions).+*/
60   |   FIELD_TENSOR_SYMMETRIC = 0x38,
61   |   /*+Shear tensor field, both symmetric and inferring last diagonal from the
62   |     opposite of the sum of the others.+*/
63   |   FIELD_TENSOR_SHEAR     = 0x39
64   | } field_plot_type;
65   | 
66   | /* Core stuff in illuminator.c */
67   | int IllDrawTet (PetscScalar *coords, PetscScalar *vals, PetscScalar isoquant,
68   | 		PetscScalar *color);
69   | int IllDrawHex (PetscScalar *coords, PetscScalar *vals, PetscScalar isoquant,
70   | 		PetscScalar *color);
71   | int IllDraw3DBlock
72   | (int xd, int yd, int zd, int xs, int ys, int zs, int xm, int ym, int zm,
73   |  PetscScalar *minmax, PetscScalar *vals, int skip,
74   |  int n_quants, PetscScalar *isoquants, PetscScalar *colors);
75   | 
76   | /* Utility stuff in utility.c */
77   | int auto_scale
78   | (PetscScalar *global_array, int points, int num_fields, int display_field,
79   |  field_plot_type fieldtype, int dimensions, PetscScalar *scale);
80   | int minmax_scale
81   | (PetscScalar *global_array, int points, int num_fields, int display_field,
82   |  field_plot_type fieldtype, int dimensions, PetscScalar *minmax);
83   | void field_indices (int nfields, int ds, field_plot_type *plottypes,
84   | 		    int *indices);
85   | 
86   | /* PETSc stuff; xcut, ycut and zcut should almost always be PETSC_FALSE */
87   | int DATriangulateRange
88   | (ISurface Surf, DA theda, Vec globalX, int this, PetscScalar *minmax,
89   |  int n_quants, PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax,
90   |  int ymin,int ymax, int zmin,int zmax);
91   | int DATriangulateLocalRange
92   | (ISurface Surf, DA theda, Vec localX, int this, PetscScalar *minmax,
93   |  int n_quants, PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax,
94   |  int ymin,int ymax, int zmin,int zmax);
95   | static inline int DATriangulate
96   | (ISurface Surf, DA theda, Vec globalX, int this, PetscScalar *minmax,
97   |  int n_quants, PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut,
98   |  PetscTruth ycut, PetscTruth zcut)
99   | { return DATriangulateRange
100  |     (Surf, theda, globalX, this, minmax, n_quants, isoquants, colors,
101  |      0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); }
102  | static inline int DATriangulateLocal
103  | (ISurface Surf, DA theda, Vec localX, int this, PetscScalar *minmax,
104  |  int n_quants, PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut,
105  |  PetscTruth ycut, PetscTruth zcut)
106  | { return DATriangulateLocalRange
107  |     (Surf, theda, localX, this, minmax, n_quants, isoquants, colors,
108  |      0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); }
109  | /* Note: IllErrorHandler is deprecated, use PETSc SETERRQ instead. */
110  | int IllErrorHandler (int id, char *message);
111  | 
112  | /* Plotting functions to render data into an RGB buffer, transform #defines */
113  | int render_scale_2d (IDisplay Disp, field_plot_type fieldtype, int symmetry);
114  | int render_composition_path
115  | (IDisplay Disp, PetscScalar *comp_array, int gridpoints, int num_fields,
116  |  field_plot_type fieldtype, PetscScalar *scale,
117  |  PetscScalar red,PetscScalar green,PetscScalar blue,PetscScalar alpha);
118  | int render_rgb_local_2d
119  | (IDisplay Disp, PetscScalar *global_array, int num_fields, int display_field,
120  |  field_plot_type fieldtype, PetscScalar *scale, int nx,int ny, int xs,int ys,
121  |  int xm,int ym, int transform, IDisplay SDisp, PetscScalar dpred,
122  |  PetscScalar dpgreen, PetscScalar dpblue, PetscScalar dpalpha);
123  | int render_rgb_local_3d
124  | (IDisplay Disp, ISurface Surf, PetscScalar *eye, PetscScalar *dir,
125  |  PetscScalar *right);
126  | #define FLIP_HORIZONTAL 0x01
127  | #define FLIP_VERTICAL   0x02
128  | #define ROTATE_LEFT     0x04
129  | 
130  | /* IlluMulti load/save stuff, including compression #defines */
131  | int IlluMultiLoad
132  | (MPI_Comm comm, char *basename, DA *theda, PetscScalar *wx,PetscScalar *wy,
133  |  PetscScalar *wz, field_plot_type **fieldtypes, int *usermetacount,
134  |  char ***usermetanames, char ***usermetadata);
135  | int IlluMultiRead
136  | (MPI_Comm comm, DA theda, Vec X, char *basename, int *usermetacount,
137  |  char ***usermetanames, char ***usermetadata);
138  | int IlluMultiSave
139  | (MPI_Comm comm, DA theda, Vec X, char *basename, PetscScalar wx,PetscScalar wy,
140  |  PetscScalar wz, field_plot_type *fieldtypes, int usermetacount,
141  |  char **usermetanames, char **usermetadata, int compressed);
142  | #define COMPRESS_INT_MASK  0x30
143  | #define COMPRESS_INT_NONE  0x00
144  | #define COMPRESS_INT_LONG  0x10
145  | #define COMPRESS_INT_SHORT 0x20
146  | #define COMPRESS_INT_CHAR  0x30
147  | #define COMPRESS_GZIP_MASK 0x0F
148  | #define COMPRESS_GZIP_NONE 0x00
149  | #define COMPRESS_GZIP_FAST 0x01
150  | #define COMPRESS_GZIP_BEST 0x0A
151  | 
152  | /* Alternate storage format import/export */
153  | typedef enum { DEFAULT=0, X_UP, X_DOWN, Y_UP, Y_DOWN, Z_UP, Z_DOWN } IllLayerStyle;
154  | typedef enum { PPM=0, TIFF, PNG, JPEG, GIF } IllImageFormat;
155  | int IllImageSave
156  | (MPI_Comm comm, DA theda, Vec X, char *basename,
157  |  int redfield,int greenfield,int bluefield,
158  |  PetscScalar *rgbmin,PetscScalar *rgbmax, int *coordrange,
159  |  IllLayerStyle layer, IllImageFormat format);
160  | 
161  | /* Geomview stuff */
162  | int GeomviewBegin (MPI_Comm comm, IDisplay *newdisp);
163  | int GeomviewEnd (MPI_Comm comm, IDisplay oldisp);
164  | int GeomviewDisplayTriangulation
165  | (MPI_Comm comm, ISurface Surf, IDisplay Disp, PetscScalar *minmax, char *name,
166  |  PetscTruth transparent);
167  | 
168  | /* Imlib2 stuff */
169  | #ifdef IMLIB2_EXISTS
170  | #include <X11/Xlib.h>
171  | #include <Imlib2.h>
172  | int imlib2_render_triangles (DATA32 *data, int width, int height,
173  | 			     int num_triangles, int *triangle_coords,
174  | 			     PetscScalar *triangle_colors, int color_skip,
175  | 			     PetscScalar *triangle_shades, int shade_skip);
176  | #endif
177  | 
178  | #endif /* ILLUMINATOR_H */