1 | /***************************************
2 | $Header: /cvsroot/petscgraphics/illuminator.h,v 1.22 2004/08/17 15:05:42 hazelsct Exp $
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 |
12 | /*+ A value of
13 | +latex+{\tt field\_plot\_type}
14 | +html+ <tt>field_plot_type</tt>
15 | is attached to each field in a simulation in order to visualize them
16 | properly. Types are as follows:
17 | +*/
18 | typedef enum {
19 | /*+Scalar field.+*/
20 | FIELD_SCALAR = 0x00,
21 | /*+Ternary composition field with two components (third component is inferred
22 | from first two).+*/
23 | FIELD_TERNARY = 0x10,
24 | /*+Vector field.+*/
25 | FIELD_VECTOR = 0x20,
26 | /*+Full ds*ds tensor field, e.g. transformation.+*/
27 | FIELD_TENSOR_FULL = 0x30,
28 | /*+Symmetric tensor field (using lines in principal stress directions).+*/
29 | FIELD_TENSOR_SYMMETRIC = 0x38,
30 | /*+Symmetric tensor field, inferring last diagonal from the opposite of the
31 | sum of the others.+*/
32 | FIELD_TENSOR_SYMMETRIC_ZERODIAG = 0x39
33 | } field_plot_type;
34 |
35 | /* Core stuff in illuminator.c */
36 | int IllDrawTet (PetscScalar *coords, PetscScalar *vals, PetscScalar isoquant,
37 | PetscScalar *color);
38 | int IllDrawHex (PetscScalar *coords, PetscScalar *vals, PetscScalar isoquant,
39 | PetscScalar *color);
40 | int IllDraw3DBlock
41 | (int xd, int yd, int zd, int xs, int ys, int zs, int xm, int ym, int zm,
42 | PetscScalar *minmax, PetscScalar *vals, int skip,
43 | int n_quants, PetscScalar *isoquants, PetscScalar *colors);
44 |
45 | /* Utility stuff in utility.c */
46 | int minmax_scale
47 | (PetscScalar *global_array, int points, int num_fields, int display_field,
48 | field_plot_type fieldtype, int dimensions, PetscScalar *minmax);
49 | void field_indices (int nfields, int ds, field_plot_type *plottypes,
50 | int *indices);
51 |
52 | /* PETSc stuff; xcut, ycut and zcut should almost always be PETSC_FALSE */
53 | int DATriangulateRange
54 | (DA theda, Vec globalX, int this, PetscScalar *minmax, int n_quants,
55 | PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin,
56 | int ymax, int zmin,int zmax);
57 | int DATriangulateLocalRange
58 | (DA theda, Vec localX, int this, PetscScalar *minmax, int n_quants,
59 | PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin,
60 | int ymax, int zmin,int zmax);
61 | static inline int DATriangulate
62 | (DA theda, Vec globalX, int this, PetscScalar *minmax, int n_quants,
63 | PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut, PetscTruth ycut,
64 | PetscTruth zcut)
65 | { DATriangulateRange (theda, globalX, this, minmax, n_quants, isoquants,
66 | colors, 0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); }
67 | static inline int DATriangulateLocal
68 | (DA theda, Vec localX, int this, PetscScalar *minmax, int n_quants,
69 | PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut, PetscTruth ycut,
70 | PetscTruth zcut)
71 | { DATriangulateLocalRange (theda, localX, this, minmax, n_quants, isoquants,
72 | colors, 0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); }
73 | int IllErrorHandler (int id, char *message);
74 |
75 | /* Plotting functions to render data into an RGB buffer */
76 | int render_rgb_local_2d
77 | (guchar *rgb, int rwidth, int rheight, int bytes_per_pixel,
78 | PetscScalar *global_array, int num_fields, int display_field,
79 | field_plot_type fieldtype, PetscScalar *minmax, int nx,int ny, int xs,int ys,
80 | int xm,int ym);
81 | int render_rgb_local_3d
82 | (guchar *rgb, int rwidth, int rheight, int bytes_per_pixel,
83 | int num_triangles, PetscScalar *vertices,
84 | PetscScalar *eye, PetscScalar *dir, PetscScalar *right);
85 |
86 | /* IlluMulti load/save stuff, including compression #defines */
87 | int IlluMultiLoad
88 | (char *basename, DA *theda, PetscScalar *wx,PetscScalar *wy,PetscScalar *wz,
89 | field_plot_type **fieldtypes, int *usermetacount, char ***usermetanames,
90 | char ***usermetadata);
91 | int IlluMultiRead (DA theda, Vec X, char *basename, int *usermetacount,
92 | char ***usermetanames, char ***usermetadata);
93 | int IlluMultiSave
94 | (DA theda, Vec X, char *basename, PetscScalar wx,PetscScalar wy,PetscScalar wz,
95 | field_plot_type *fieldtypes, int usermetacount, char **usermetanames,
96 | char **usermetadata, int compressed);
97 | #define COMPRESS_INT_MASK 0x30
98 | #define COMPRESS_INT_NONE 0x00
99 | #define COMPRESS_INT_LONG 0x10
100 | #define COMPRESS_INT_SHORT 0x20
101 | #define COMPRESS_INT_CHAR 0x30
102 | #define COMPRESS_GZIP_MASK 0x0F
103 | #define COMPRESS_GZIP_NONE 0x00
104 | #define COMPRESS_GZIP_FAST 0x01
105 | #define COMPRESS_GZIP_BEST 0x0A
106 |
107 | /* Geomview stuff */
108 | int GeomviewBegin (MPI_Comm comm);
109 | int GeomviewEnd (MPI_Comm comm);
110 | int GeomviewDisplayTriangulation (MPI_Comm comm, PetscScalar *minmax,
111 | char *name, PetscTruth transparent);
112 |
113 | #endif /* ILLUMINATOR_H */