1 | /*************************************** 2 | $Header: /cvsroot/petscgraphics/illuminator.h,v 1.27 2004/11/04 22:47:45 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 | /*+Shear tensor field, both symmetric and inferring last diagonal from the 31 | opposite of the sum of the others.+*/ 32 | FIELD_TENSOR_SHEAR = 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 auto_scale 47 | (PetscScalar *global_array, int points, int num_fields, int display_field, 48 | field_plot_type fieldtype, int dimensions, PetscScalar *scale); 49 | int minmax_scale 50 | (PetscScalar *global_array, int points, int num_fields, int display_field, 51 | field_plot_type fieldtype, int dimensions, PetscScalar *minmax); 52 | void field_indices (int nfields, int ds, field_plot_type *plottypes, 53 | int *indices); 54 | 55 | /* PETSc stuff; xcut, ycut and zcut should almost always be PETSC_FALSE */ 56 | int DATriangulateRange 57 | (DA theda, Vec globalX, int this, PetscScalar *minmax, int n_quants, 58 | PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin, 59 | int ymax, int zmin,int zmax); 60 | int DATriangulateLocalRange 61 | (DA theda, Vec localX, int this, PetscScalar *minmax, int n_quants, 62 | PetscScalar *isoquants, PetscScalar *colors, int xmin,int xmax, int ymin, 63 | int ymax, int zmin,int zmax); 64 | static inline int DATriangulate 65 | (DA theda, Vec globalX, int this, PetscScalar *minmax, int n_quants, 66 | PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut, PetscTruth ycut, 67 | PetscTruth zcut) 68 | { DATriangulateRange (theda, globalX, this, minmax, n_quants, isoquants, 69 | colors, 0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); } 70 | static inline int DATriangulateLocal 71 | (DA theda, Vec localX, int this, PetscScalar *minmax, int n_quants, 72 | PetscScalar *isoquants, PetscScalar *colors, PetscTruth xcut, PetscTruth ycut, 73 | PetscTruth zcut) 74 | { DATriangulateLocalRange (theda, localX, this, minmax, n_quants, isoquants, 75 | colors, 0,xcut?-2:-1, 0,ycut?-2:-1, 0,zcut?-2:-1); } 76 | int IllErrorHandler (int id, char *message); 77 | 78 | /* Plotting functions to render data into an RGB buffer */ 79 | int render_scale_2d (guchar *rgb, int rwidth, int rheight, int bytes_per_pixel, 80 | field_plot_type fieldtype, int symmetry); 81 | int render_rgb_local_2d 82 | (guchar *rgb, int rwidth, int rheight, int rowskip, int bytes_per_pixel, 83 | PetscScalar *global_array, int num_fields, int display_field, 84 | field_plot_type fieldtype, PetscScalar *minmax, int nx,int ny, int xs,int ys, 85 | int xm,int ym); 86 | int render_rgb_local_3d 87 | (guchar *rgb, int rwidth, int rheight, int bytes_per_pixel, 88 | int num_triangles, PetscScalar *vertices, 89 | PetscScalar *eye, PetscScalar *dir, PetscScalar *right); 90 | 91 | /* IlluMulti load/save stuff, including compression #defines */ 92 | int IlluMultiLoad 93 | (char *basename, DA *theda, PetscScalar *wx,PetscScalar *wy,PetscScalar *wz, 94 | field_plot_type **fieldtypes, int *usermetacount, char ***usermetanames, 95 | char ***usermetadata); 96 | int IlluMultiRead (DA theda, Vec X, char *basename, int *usermetacount, 97 | char ***usermetanames, char ***usermetadata); 98 | int IlluMultiSave 99 | (DA theda, Vec X, char *basename, PetscScalar wx,PetscScalar wy,PetscScalar wz, 100 | field_plot_type *fieldtypes, int usermetacount, char **usermetanames, 101 | char **usermetadata, int compressed); 102 | #define COMPRESS_INT_MASK 0x30 103 | #define COMPRESS_INT_NONE 0x00 104 | #define COMPRESS_INT_LONG 0x10 105 | #define COMPRESS_INT_SHORT 0x20 106 | #define COMPRESS_INT_CHAR 0x30 107 | #define COMPRESS_GZIP_MASK 0x0F 108 | #define COMPRESS_GZIP_NONE 0x00 109 | #define COMPRESS_GZIP_FAST 0x01 110 | #define COMPRESS_GZIP_BEST 0x0A 111 | 112 | /* Geomview stuff */ 113 | int GeomviewBegin (MPI_Comm comm); 114 | int GeomviewEnd (MPI_Comm comm); 115 | int GeomviewDisplayTriangulation (MPI_Comm comm, PetscScalar *minmax, 116 | char *name, PetscTruth transparent); 117 | 118 | /* Imlib2 stuff */ 119 | #ifdef IMLIB2_EXISTS 120 | #include <X11/Xlib.h> 121 | #include <Imlib2.h> 122 | int imlib2_render_triangles (DATA32 *data, int width, int height, 123 | int num_triangles, int *triangle_coords, 124 | PetscScalar *triangle_colors, int color_skip, 125 | PetscScalar *triangle_shades, int shade_skip); 126 | #endif 127 | 128 | #endif /* ILLUMINATOR_H */