Salome HOME
IMP 0017328: min and max scalar map of results given at gauss points. A fix for the...
[modules/visu.git] / src / PIPELINE / VISU_OpenGLPointSpriteMapper.cxx
1 //  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 //  This library is free software; you can redistribute it and/or
7 //  modify it under the terms of the GNU Lesser General Public
8 //  License as published by the Free Software Foundation; either
9 //  version 2.1 of the License.
10 //
11 //  This library is distributed in the hope that it will be useful,
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 //  Lesser General Public License for more details.
15 //
16 //  You should have received a copy of the GNU Lesser General Public
17 //  License along with this library; if not, write to the Free Software
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //  VISU OBJECT : interactive object for VISU entities implementation
23 // File:    VISU_OpenGLPointSpriteMapper.cxx
24 // Author:  
25 // Module : VISU
26 //
27 #include "VISU_OpenGLPointSpriteMapper.hxx"
28
29 #include "SVTK_Extension.h"
30
31 #include <vtkCamera.h>
32 #include <vtkCellArray.h>
33 #include <vtkCellData.h>
34 #include <vtkCommand.h>
35 #include <vtkImageData.h>
36 #include <vtkMatrix4x4.h>
37 #include <vtkObjectFactory.h>
38 #include <vtkRenderer.h>
39 #include <vtkRenderWindow.h>
40 #include <vtkPointData.h>
41 #include <vtkPolyData.h>
42 #include <vtkPolygon.h>
43 #include <vtkProperty.h>
44 #include <vtkTimerLog.h>
45 #include <vtkTriangle.h>
46
47 #include <stdio.h>
48 #include <cmath>
49 #include <string>
50
51 #include "utilities.h"
52
53 #ifndef WNT
54 # ifndef GLX_GLXEXT_LEGACY
55 #  define GLX_GLXEXT_LEGACY
56 # endif
57 # include <GL/glx.h>
58 # include <dlfcn.h>
59 #else
60 # include <wingdi.h>
61 #endif
62
63 #ifndef VTK_IMPLEMENT_MESA_CXX
64 vtkCxxRevisionMacro(VISU_OpenGLPointSpriteMapper, "Revision$");
65 vtkStandardNewMacro(VISU_OpenGLPointSpriteMapper);
66 #endif
67
68 // some definitions for what the polydata has in it
69 #define VTK_PDPSM_COLORS             0x0001
70 #define VTK_PDPSM_CELL_COLORS        0x0002
71 #define VTK_PDPSM_POINT_TYPE_FLOAT   0x0004
72 #define VTK_PDPSM_POINT_TYPE_DOUBLE  0x0008
73 #define VTK_PDPSM_NORMAL_TYPE_FLOAT  0x0010
74 #define VTK_PDPSM_NORMAL_TYPE_DOUBLE 0x0020
75 #define VTK_PDPSM_OPAQUE_COLORS      0x0040
76
77 #ifndef APIENTRY
78 #define APIENTRY
79 #endif
80 #ifndef APIENTRYP
81 #define APIENTRYP APIENTRY *
82 #endif
83
84 #ifndef GL_OBJECT_INFO_LOG_LENGTH_ARB
85 #define GL_OBJECT_INFO_LOG_LENGTH_ARB     0x8B84
86 #endif
87
88 #ifndef GL_VERTEX_SHADER_ARB
89 #define GL_VERTEX_SHADER_ARB              0x8B31
90 #endif
91
92 #ifndef GL_VERTEX_PROGRAM_POINT_SIZE_ARB
93 #define GL_VERTEX_PROGRAM_POINT_SIZE_ARB  0x8642
94 #endif
95
96 #ifndef GL_ARB_point_sprite
97 #define GL_POINT_SPRITE_ARB               0x8861
98 #define GL_COORD_REPLACE_ARB              0x8862
99 #endif
100
101 #ifndef GL_ARB_shader_objects
102 typedef char GLcharARB;
103 #endif
104
105 #ifndef GL_ARB_vertex_buffer_object
106 typedef ptrdiff_t GLsizeiptrARB;
107
108 #define GL_ARRAY_BUFFER_ARB               0x8892
109 #define GL_STATIC_DRAW_ARB                0x88E4
110 #endif
111
112 typedef void (APIENTRYP PFNGLSHADERSOURCEARBPROC) (GLhandleARB shaderObj, GLsizei count, const GLcharARB* *string, const GLint *length);
113 typedef GLhandleARB (APIENTRYP PFNGLCREATESHADEROBJECTARBPROC) (GLenum shaderType);
114 typedef void (APIENTRYP PFNGLCOMPILESHADERARBPROC) (GLhandleARB shaderObj);
115 typedef GLhandleARB (APIENTRYP PFNGLCREATEPROGRAMOBJECTARBPROC) (void);
116 typedef void (APIENTRYP PFNGLATTACHOBJECTARBPROC) (GLhandleARB containerObj, GLhandleARB obj);
117 typedef void (APIENTRYP PFNGLLINKPROGRAMARBPROC) (GLhandleARB programObj);
118 typedef void (APIENTRYP PFNGLUSEPROGRAMOBJECTARBPROC) (GLhandleARB programObj);
119 typedef void (APIENTRYP PFNGLGETOBJECTPARAMETERIVARBPROC) (GLhandleARB obj, GLenum pname, GLint *params);
120 typedef void (APIENTRYP PFNGLGETINFOLOGARBPROC) (GLhandleARB obj, GLsizei maxLength, GLsizei *length, GLcharARB *infoLog);
121 typedef GLint (APIENTRYP PFNGLGETATTRIBLOCATIONARBPROC) (GLhandleARB programObj, const GLcharARB *name);
122 typedef void (APIENTRYP PFNGLVERTEXATTRIB1FARBPROC) (GLuint index, GLfloat x);
123
124 typedef void (APIENTRYP PFNGLBINDBUFFERARBPROC) (GLenum target, GLuint buffer);
125 typedef void (APIENTRYP PFNGLDELETEBUFFERSARBPROC) (GLsizei n, const GLuint *buffers);
126 typedef void (APIENTRYP PFNGLGENBUFFERSARBPROC) (GLsizei n, GLuint *buffers);
127 typedef void (APIENTRYP PFNGLBUFFERDATAARBPROC) (GLenum target, GLsizeiptrARB size, const GLvoid *data, GLenum usage);
128
129 static PFNGLSHADERSOURCEARBPROC             vglShaderSourceARB            = NULL;
130 static PFNGLCREATESHADEROBJECTARBPROC       vglCreateShaderObjectARB      = NULL;
131 static PFNGLCOMPILESHADERARBPROC            vglCompileShaderARB           = NULL;
132 static PFNGLCREATEPROGRAMOBJECTARBPROC      vglCreateProgramObjectARB     = NULL;
133 static PFNGLATTACHOBJECTARBPROC             vglAttachObjectARB            = NULL;
134 static PFNGLLINKPROGRAMARBPROC              vglLinkProgramARB             = NULL;
135 static PFNGLUSEPROGRAMOBJECTARBPROC         vglUseProgramObjectARB        = NULL;
136 static PFNGLGETOBJECTPARAMETERIVARBPROC     vglGetObjectParameterivARB    = NULL;
137 static PFNGLGETINFOLOGARBPROC               vglGetInfoLogARB              = NULL;
138 static PFNGLGETATTRIBLOCATIONARBPROC        vglGetAttribLocationARB       = NULL;
139 static PFNGLVERTEXATTRIB1FARBPROC           vglVertexAttrib1fARB          = NULL;
140
141 static PFNGLGENBUFFERSARBPROC               vglGenBuffersARB              = NULL;
142 static PFNGLBINDBUFFERARBPROC               vglBindBufferARB              = NULL;
143 static PFNGLBUFFERDATAARBPROC               vglBufferDataARB              = NULL;
144 static PFNGLDELETEBUFFERSARBPROC            vglDeleteBuffersARB           = NULL;
145
146 #ifndef WNT
147 #define GL_GetProcAddress( x )   glXGetProcAddressARB( (const GLubyte*)x )
148 #else
149 #define GL_GetProcAddress( x )   wglGetProcAddress( (const LPCSTR)x )
150 #endif
151
152 bool InitializeARB()
153 {
154   vglShaderSourceARB = (PFNGLSHADERSOURCEARBPROC)GL_GetProcAddress( "glShaderSourceARB" );
155   if( !vglShaderSourceARB )
156     return false;
157
158   vglCreateShaderObjectARB = (PFNGLCREATESHADEROBJECTARBPROC)GL_GetProcAddress( "glCreateShaderObjectARB" );
159   if( !vglCreateShaderObjectARB )
160     return false;
161
162   vglCompileShaderARB = (PFNGLCOMPILESHADERARBPROC)GL_GetProcAddress( "glCompileShaderARB" );
163   if( !vglCompileShaderARB )
164     return false;
165
166   vglCreateProgramObjectARB = (PFNGLCREATEPROGRAMOBJECTARBPROC)GL_GetProcAddress( "glCreateProgramObjectARB" );
167   if( !vglCreateProgramObjectARB )
168     return false;
169
170   vglAttachObjectARB = (PFNGLATTACHOBJECTARBPROC)GL_GetProcAddress( "glAttachObjectARB" );
171   if( !vglAttachObjectARB )
172     return false;
173
174   vglLinkProgramARB = (PFNGLLINKPROGRAMARBPROC)GL_GetProcAddress( "glLinkProgramARB" );
175   if( !vglLinkProgramARB )
176     return false;
177
178   vglUseProgramObjectARB = (PFNGLUSEPROGRAMOBJECTARBPROC)GL_GetProcAddress( "glUseProgramObjectARB" );
179   if( !vglUseProgramObjectARB )
180     return false;
181
182   vglGetObjectParameterivARB = (PFNGLGETOBJECTPARAMETERIVARBPROC)GL_GetProcAddress( "glGetObjectParameterivARB" );
183   if( !vglGetObjectParameterivARB )
184     return false;
185
186   vglGetInfoLogARB = (PFNGLGETINFOLOGARBPROC)GL_GetProcAddress( "glGetInfoLogARB" );
187   if( !vglGetInfoLogARB )
188     return false;
189
190   vglGetAttribLocationARB = (PFNGLGETATTRIBLOCATIONARBPROC)GL_GetProcAddress( "glGetAttribLocationARB" );
191   if( !vglGetAttribLocationARB )
192     return false;
193
194   vglVertexAttrib1fARB = (PFNGLVERTEXATTRIB1FARBPROC)GL_GetProcAddress( "glVertexAttrib1fARB" );
195   if( !vglVertexAttrib1fARB )
196     return false;
197
198   vglGenBuffersARB = (PFNGLGENBUFFERSARBPROC)GL_GetProcAddress( "glGenBuffersARB" );
199   if( !vglGenBuffersARB )
200     return false;
201
202   vglBindBufferARB = (PFNGLBINDBUFFERARBPROC)GL_GetProcAddress( "glBindBufferARB" );
203   if( !vglBindBufferARB )
204     return false;
205
206   vglBufferDataARB = (PFNGLBUFFERDATAARBPROC)GL_GetProcAddress( "glBufferDataARB" );
207   if( !vglBufferDataARB )
208     return false;
209
210   vglDeleteBuffersARB = (PFNGLDELETEBUFFERSARBPROC)GL_GetProcAddress( "glDeleteBuffersARB" );
211   if( !vglDeleteBuffersARB )
212     return false;
213
214   return true;
215 };
216
217 static bool IsARBInitialized = InitializeARB();
218 static vtkFloatingPointType Tolerance = 1.0 / VTK_LARGE_FLOAT;
219
220 //-----------------------------------------------------------------------------
221 // Construct empty object.
222 VISU_OpenGLPointSpriteMapper::VISU_OpenGLPointSpriteMapper()
223 {
224   this->RenderMode               = VISU_OpenGLPointSpriteMapper::Occlude;
225
226   this->ListId                   = 0;
227   this->TotalCells               = 0;
228   this->ExtensionsInitialized    = 0;
229   this->DefaultPointSize         = 20.0;
230   this->AverageCellSize          = 0.0;
231
232   this->UsePointSprites          = true;
233   this->UseTextures              = true;
234   this->UseShader                = true;
235
236   this->PrimitiveType            = VISU_OpenGLPointSpriteMapper::PointSprite;
237
238   this->PointSpriteMode          = 0;
239
240   this->PointSpriteClamp         = 256.0;
241   this->PointSpriteSize          = 0.2;
242   this->PointSpriteMinSize       = 0.1;
243   this->PointSpriteMaxSize       = 0.3;
244   this->PointSpriteMagnification = 1.0;
245
246   this->PointSpriteAlphaThreshold = 0.5;
247   this->PointSpriteOpacity       = 1.0;
248   this->PointSpriteTexture       = 0;
249
250   this->UseOpenGLMapper          = false;
251 }
252 //-----------------------------------------------------------------------------
253 VISU_OpenGLPointSpriteMapper::~VISU_OpenGLPointSpriteMapper()
254 {
255   if( PointSpriteTexture>0 )
256     glDeleteTextures( 1, &PointSpriteTexture );
257
258   if( this->LastWindow )
259     this->ReleaseGraphicsResources(this->LastWindow);
260 }
261
262 //-----------------------------------------------------------------------------
263 char* readFromFile( std::string fileName )
264 {
265   FILE* file = fopen( fileName.c_str(), "r" );
266
267   char* content = NULL;
268   int count = 0;
269
270   if( file != NULL )
271   {
272     fseek( file, 0, SEEK_END );
273     count = ftell( file );
274     rewind( file );
275
276     if( count > 0 )
277     {
278       content = ( char* )malloc( sizeof( char ) * ( count + 1 ) );
279       count = fread( content, sizeof( char ), count, file );
280       content[ count ] = '\0';
281     }
282     fclose( file );
283   }
284
285   return content;
286 }
287 //-----------------------------------------------------------------------------
288 void VISU_OpenGLPointSpriteMapper::PrintInfoLog( GLhandleARB obj )
289 {
290   int infologLength = 0;
291   int charsWritten  = 0;
292   char* infoLog;
293
294   vglGetObjectParameterivARB( obj, GL_OBJECT_INFO_LOG_LENGTH_ARB, &infologLength );
295
296   if( infologLength > 0 )
297   {
298     infoLog = ( char* )malloc( infologLength );
299     vglGetInfoLogARB( obj, infologLength, &charsWritten, infoLog );
300     printf( "%s\n", infoLog );
301     free( infoLog );
302   }
303 }
304 //-----------------------------------------------------------------------------
305 void VISU_OpenGLPointSpriteMapper::InitShader()
306 {
307   //cout << "Initializing vertex program" << endl;
308
309   std::string fileName = std::string( getenv( "VISU_ROOT_DIR") ) +
310                          "/share/salome/resources/visu/Vertex_Program_ARB.txt";
311
312   char* shader = readFromFile( fileName );
313
314   GLhandleARB VertexShader = vglCreateShaderObjectARB( GL_VERTEX_SHADER_ARB );
315   vglShaderSourceARB( VertexShader, 1, (const GLcharARB**)&shader, NULL );
316   vglCompileShaderARB( VertexShader );
317   //this->PrintInfoLog( VertexShader );
318
319   this->VertexProgram = vglCreateProgramObjectARB();
320   vglAttachObjectARB( this->VertexProgram, VertexShader );
321
322   vglLinkProgramARB( this->VertexProgram );
323   //this->PrintInfoLog( VertexProgram );
324   /*
325   cout << "Shader from " << fileName << endl;
326   for( int i = 0; i < strlen( shader ); i++ )
327     cout << shader[i];
328   cout << endl;
329
330   if( glGetError() == GL_NO_ERROR )
331     cout << "Loading vertex program... ok" << endl << endl;
332   else
333     cout << "Loading vertex program... failed" << endl << endl;
334   */
335   free( shader );
336 }
337
338 //-----------------------------------------------------------------------------
339 void VISU_OpenGLPointSpriteMapper::SetShaderVariable( const char* variable, float value )
340 {
341   //cout << this->VertexProgram << " ";
342   //cout << vglGetAttribLocationARB( this->VertexProgram, variable ) << " ";
343   //cout << variable << " " << value << endl;
344
345   vglVertexAttrib1fARB( vglGetAttribLocationARB( this->VertexProgram, variable ), value );
346 }
347
348 //-----------------------------------------------------------------------------
349 void VISU_OpenGLPointSpriteMapper::SetPrimitiveType( int thePrimitiveType )
350 {
351   if( this->PrimitiveType == thePrimitiveType )
352     return;
353
354   this->PrimitiveType = thePrimitiveType;
355   this->Modified();
356 }
357
358 //-----------------------------------------------------------------------------
359 void VISU_OpenGLPointSpriteMapper::SetPointSpriteMode( int theMode )
360 {
361   if( this->PointSpriteMode == theMode )
362     return;
363
364   this->PointSpriteMode = theMode;
365   this->Modified();
366 }
367
368 //-----------------------------------------------------------------------------
369 void VISU_OpenGLPointSpriteMapper::SetPointSpriteClamp( float theClamp )
370 {
371   if( fabs( this->PointSpriteClamp - theClamp ) < Tolerance )
372     return;
373
374   this->PointSpriteClamp = theClamp;
375   this->Modified();
376 }
377
378 //-----------------------------------------------------------------------------
379 void VISU_OpenGLPointSpriteMapper::SetAverageCellSize(float theSize)
380 {
381   if( fabs( this->AverageCellSize - theSize ) < Tolerance )
382     return;
383
384   this->AverageCellSize = theSize;
385   this->Modified();
386 }
387
388 //-----------------------------------------------------------------------------
389 void VISU_OpenGLPointSpriteMapper::SetPointSpriteSize( float theSize )
390 {
391   if( fabs( this->PointSpriteSize - theSize ) < Tolerance )
392     return;
393
394   this->PointSpriteSize = theSize;
395   this->Modified();
396 }
397
398 //-----------------------------------------------------------------------------
399 void VISU_OpenGLPointSpriteMapper::SetPointSpriteMinSize( float theMinSize )
400 {
401   if( fabs( this->PointSpriteMinSize - theMinSize ) < Tolerance )
402     return;
403
404   this->PointSpriteMinSize = theMinSize;
405   this->Modified();
406 }
407
408 //-----------------------------------------------------------------------------
409 void VISU_OpenGLPointSpriteMapper::SetPointSpriteMaxSize( float theMaxSize )
410 {
411   if( fabs( this->PointSpriteMaxSize - theMaxSize ) < Tolerance )
412     return;
413
414   this->PointSpriteMaxSize = theMaxSize;
415   this->Modified();
416 }
417
418 //-----------------------------------------------------------------------------
419 void VISU_OpenGLPointSpriteMapper::SetPointSpriteMagnification( float theMagnification )
420 {
421   if( fabs( this->PointSpriteMagnification - theMagnification ) < Tolerance )
422     return;
423
424   this->PointSpriteMagnification = theMagnification;
425   this->Modified();
426 }
427
428 //-----------------------------------------------------------------------------
429 void VISU_OpenGLPointSpriteMapper::SetPointSpriteAlphaThreshold( float theAlphaThreshold )
430 {
431   if( fabs( this->PointSpriteAlphaThreshold - theAlphaThreshold ) < Tolerance )
432     return;
433
434   this->PointSpriteAlphaThreshold = theAlphaThreshold;
435   this->Modified();
436 }
437
438 //-----------------------------------------------------------------------------
439 bool VISU_OpenGLPointSpriteMapper::InitExtensions()
440 {
441   if( this->ExtensionsInitialized )
442     return true;
443
444   InitializeARB();
445
446   char* ext = (char*)glGetString( GL_EXTENSIONS );
447   //cout << "OpenGL extensions : " << ext << endl;
448
449   if( !IsARBInitialized ||
450       strstr( ext, "GL_ARB_point_sprite" ) == NULL ||
451       strstr( ext, "GL_ARB_shader_objects" ) == NULL ||
452       strstr( ext, "GL_ARB_vertex_buffer_object" ) == NULL )
453   {
454     vtkWarningMacro(<<"Initializing ARB extensions failed");
455     this->UseOpenGLMapper = true;
456
457     return false;
458   }
459
460   if( this->UseShader )
461     this->InitShader();
462
463   this->ExtensionsInitialized = 1;
464   return true;
465 }
466
467 //-----------------------------------------------------------------------------
468 float ViewToDisplay( vtkRenderer* theRenderer )
469 {
470   vtkFloatingPointType p1[3], p2[3];
471
472   theRenderer->SetViewPoint( 0.0, 0.0, 0.0 );
473   theRenderer->ViewToDisplay();
474   theRenderer->GetDisplayPoint( p1 );
475
476   theRenderer->SetViewPoint( 1.0, 1.0, 1.0 );
477   theRenderer->ViewToDisplay();
478   theRenderer->GetDisplayPoint( p2 );
479
480   vtkFloatingPointType coefficient = sqrt( pow( p2[0] - p1[0], 2 ) + pow( p2[1] - p1[1], 2 ) ) / sqrt( 2. );
481   //cout << p1[0] << " " << p1[1] << " " << p1[2] << endl;
482   //cout << p2[0] << " " << p2[1] << " " << p2[2] << endl;
483   //cout << "ZOOM  : " << coefficient << endl;
484
485   return coefficient;
486 }
487
488 //-----------------------------------------------------------------------------
489 //
490 // Receives from Actor -> maps data to primitives
491 //
492 void VISU_OpenGLPointSpriteMapper::RenderPiece(vtkRenderer *ren, vtkActor *act)
493 {
494   bool isUseThisMapper = this->PrimitiveType != VISU_OpenGLPointSpriteMapper::GeomSphere;
495
496   if( isUseThisMapper )
497     this->InitExtensions();
498
499   if( !isUseThisMapper )
500   {
501     MAPPER_SUPERCLASS::RenderPiece( ren, act );
502     return;
503   }
504
505   vtkIdType numPts;
506   vtkPolyData *input= this->GetInput();
507
508   //
509   // make sure that we've been properly initialized
510   //
511   if (ren->GetRenderWindow()->CheckAbortStatus())
512     return;
513
514   if ( input == NULL )
515   {
516     vtkErrorMacro(<< "No input!");
517     return;
518   }
519   else
520   {
521     this->InvokeEvent(vtkCommand::StartEvent,NULL);
522     input->Update();
523     this->InvokeEvent(vtkCommand::EndEvent,NULL);
524
525     numPts = input->GetNumberOfPoints();
526   }
527
528   if (numPts == 0)
529   {
530     vtkDebugMacro(<< "No points!");
531     return;
532   }
533
534   if ( this->LookupTable == NULL )
535     this->CreateDefaultLookupTable();
536
537   // make sure our window is current
538   ren->GetRenderWindow()->MakeCurrent();
539
540   if( this->UsePointSprites ) //&& this->PrimitiveType == VISU_OpenGLPointSpriteMapper::PointSprite )
541     this->InitPointSprites();
542
543   // Initializing the texture for Point Sprites
544   if( this->UseTextures && this->PrimitiveType == VISU_OpenGLPointSpriteMapper::PointSprite )
545     this->InitTextures();
546
547   vglUseProgramObjectARB( this->VertexProgram );
548   float aViewToDisplay = ViewToDisplay( ren );
549   this->SetShaderVariable( "mode",          this->PointSpriteMode );
550   this->SetShaderVariable( "clampSize",     this->PointSpriteClamp );
551   this->SetShaderVariable( "geomSize",      aViewToDisplay * this->AverageCellSize * this->PointSpriteSize );
552   this->SetShaderVariable( "minSize",       aViewToDisplay * this->AverageCellSize * this->PointSpriteMinSize );
553   this->SetShaderVariable( "maxSize",       aViewToDisplay * this->AverageCellSize * this->PointSpriteMaxSize );
554   this->SetShaderVariable( "magnification", this->PointSpriteMagnification );
555
556   //
557   // if something has changed regenerate colors and display lists
558   // if required
559   //
560   int noAbort=1;
561   if ( this->GetMTime() > this->BuildTime ||
562        input->GetMTime() > this->BuildTime ||
563        act->GetProperty()->GetMTime() > this->BuildTime ||
564        ren->GetRenderWindow() != this->LastWindow)
565   {
566 #ifdef _DEBUG_RENDERING_PERFORMANCE_
567     // To control when the mapper is recalculated
568     MESSAGE( "VISU_OpenGLPointSpriteMapper::RenderPiece - "
569              <<(this->GetMTime() > this->BuildTime)<<"; "
570              <<(input->GetMTime() > this->BuildTime)<<"; "
571              <<(act->GetProperty()->GetMTime() > this->BuildTime)<<"; ");
572 #endif
573     // sets this->Colors as side effect
574     this->MapScalars( act->GetProperty()->GetOpacity() );
575
576     if (!this->ImmediateModeRendering &&
577         !this->GetGlobalImmediateModeRendering())
578     {
579       this->ReleaseGraphicsResources(ren->GetRenderWindow());
580       this->LastWindow = ren->GetRenderWindow();
581
582       // get a unique display list id
583       this->ListId = glGenLists(1);
584       glNewList(this->ListId,GL_COMPILE);
585
586       noAbort = this->Draw(ren,act);
587       glEndList();
588
589       // Time the actual drawing
590       this->Timer->StartTimer();
591       glCallList(this->ListId);
592       this->Timer->StopTimer();
593     }
594     else
595     {
596       this->ReleaseGraphicsResources(ren->GetRenderWindow());
597       this->LastWindow = ren->GetRenderWindow();
598     }
599     if (noAbort)
600       this->BuildTime.Modified();
601   }
602   // if nothing changed but we are using display lists, draw it
603   else
604   {
605     if (!this->ImmediateModeRendering &&
606         !this->GetGlobalImmediateModeRendering())
607     {
608       // Time the actual drawing
609       this->Timer->StartTimer();
610       glCallList(this->ListId);
611       this->Timer->StopTimer();
612     }
613   }
614
615   // if we are in immediate mode rendering we always
616   // want to draw the primitives here
617   if (this->ImmediateModeRendering ||
618       this->GetGlobalImmediateModeRendering())
619   {
620     // sets this->Colors as side effect
621     this->MapScalars( act->GetProperty()->GetOpacity() );
622
623     // Time the actual drawing
624     this->Timer->StartTimer();
625     this->Draw(ren,act);
626     this->Timer->StopTimer();
627   }
628
629   this->TimeToDraw = (float)this->Timer->GetElapsedTime();
630
631   // If the timer is not accurate enough, set it to a small
632   // time so that it is not zero
633   if ( this->TimeToDraw == 0.0 )
634     this->TimeToDraw = 0.0001;
635
636   vglUseProgramObjectARB( 0 );
637
638   if( this->UsePointSprites ) //&& this->PrimitiveType == VISU_OpenGLPointSpriteMapper::PointSprite )
639     this->CleanupPointSprites();
640 }
641
642 //-----------------------------------------------------------------------------
643 float VISU_OpenGLPointSpriteMapper::GetMaximumSupportedSize()
644 {
645   float maximumSupportedSize = 512.0;
646   //glGetFloatv( GL_POINT_SIZE_MAX_ARB, &maximumSupportedSize );
647
648   return maximumSupportedSize;
649 }
650
651 //-----------------------------------------------------------------------------
652 void VISU_OpenGLPointSpriteMapper::InitPointSprites()
653 {
654   glEnable( GL_POINT_SPRITE_ARB );
655   glEnable( GL_VERTEX_PROGRAM_POINT_SIZE_ARB );
656
657   glPushAttrib( GL_COLOR_BUFFER_BIT | GL_CURRENT_BIT | GL_DEPTH_BUFFER_BIT | GL_ENABLE_BIT | GL_LIGHTING_BIT );
658
659   this->RenderMode = this->PointSpriteOpacity < 1.0 ? VISU_OpenGLPointSpriteMapper::Accumulate : VISU_OpenGLPointSpriteMapper::Occlude;
660
661   switch (this->RenderMode)
662   {
663     case VISU_OpenGLPointSpriteMapper::Accumulate:
664     {
665       glDepthFunc( GL_LESS );
666       glEnable( GL_DEPTH_TEST );
667
668       glEnable( GL_BLEND );
669       glBlendFunc( GL_SRC_ALPHA, GL_ONE );
670
671       glEnable( GL_ALPHA_TEST );
672       glAlphaFunc( GL_GREATER, this->PointSpriteAlphaThreshold );
673       break;
674     }
675
676     case VISU_OpenGLPointSpriteMapper::Occlude:
677     {
678       glDepthFunc( GL_LEQUAL );
679       glEnable( GL_DEPTH_TEST );
680
681       glEnable( GL_ALPHA_TEST );
682       glAlphaFunc( GL_GREATER, this->PointSpriteAlphaThreshold );
683
684       glDisable( GL_BLEND );
685       break;
686     }
687
688     default:
689     {
690       break;
691     }
692   }
693   // Disable Lighting/Shading.
694   glDisable( GL_LIGHTING );
695
696   // Disable material properties
697   glDisable( GL_COLOR_MATERIAL );
698 }
699
700 //-----------------------------------------------------------------------------
701 void VISU_OpenGLPointSpriteMapper::CleanupPointSprites()
702 {
703   // Set GL params back to normal to stop other vtkMappers displaying wrongly
704   glPopAttrib();
705
706   glDisable( GL_VERTEX_PROGRAM_POINT_SIZE_ARB );
707   glDisable( GL_POINT_SPRITE_ARB );
708 }
709
710
711 //-----------------------------------------------------------------------------
712 void
713 VISU_OpenGLPointSpriteMapper
714 ::SetImageData( vtkImageData* theImageData )
715 {
716   if(GetImageData() == theImageData)
717     return;
718   this->ImageData = theImageData;
719   this->Modified();
720 }
721
722 vtkImageData*
723 VISU_OpenGLPointSpriteMapper
724 ::GetImageData()
725 {
726   return this->ImageData.GetPointer();
727 }
728
729
730 //-----------------------------------------------------------------------------
731 void VISU_OpenGLPointSpriteMapper::InitTextures()
732 {
733   //cout << "VISU_OpenGLPointSpriteMapper::InitTextures " << this->GetImageData() << endl;
734   if( !this->GetImageData() )
735     return;
736
737   glTexParameteri( GL_TEXTURE_2D, GL_TEXTURE_WRAP_S, GL_REPEAT );
738   glTexParameteri( GL_TEXTURE_2D, GL_TEXTURE_WRAP_T, GL_REPEAT );
739   glTexParameterf( GL_TEXTURE_2D, GL_TEXTURE_MIN_FILTER, GL_LINEAR );
740   glTexParameterf( GL_TEXTURE_2D, GL_TEXTURE_MAG_FILTER, GL_LINEAR );
741
742   int* aSize = GetImageData()->GetDimensions();
743   unsigned char* dataPtr = (unsigned char*)GetImageData()->GetScalarPointer();
744   glTexImage2D( GL_TEXTURE_2D, 0, GL_RGBA, aSize[0], aSize[1], 0,
745                 GL_RGBA, GL_UNSIGNED_BYTE, dataPtr );
746
747   //glClear( GL_COLOR_BUFFER_BIT | GL_DEPTH_BUFFER_BIT );
748   glEnable( GL_TEXTURE_2D );
749   glTexEnvf( GL_POINT_SPRITE_ARB, GL_COORD_REPLACE_ARB, GL_TRUE );
750   glBindTexture( GL_TEXTURE_2D, this->PointSpriteTexture );
751 }
752
753
754 //-----------------------------------------------------------------------------
755 int ComputeHue( int r, int g, int b )
756 {
757   int h = 0;
758
759   int max = r;
760   int whatmax = 0;
761   if( g > max ) {
762     max = g;
763     whatmax = 1;
764   }
765   if( b > max ) {
766     max = b;
767     whatmax = 2;
768   }
769
770   int min = r;
771   if ( g < min ) min = g;
772   if ( b < min ) min = b;
773   int delta = max-min;
774
775   if( delta == 0 )
776     h = 0;
777   else if( whatmax == 0 ) {
778     if ( g >= b )
779       h = (120*(g-b)+delta)/(2*delta);
780     else
781       h = (120*(g-b+delta)+delta)/(2*delta) + 300;
782   }
783   else if( whatmax == 1 ) {
784     if ( b > r )
785       h = 120 + (120*(b-r)+delta)/(2*delta);
786     else
787       h = 60 + (120*(b-r+delta)+delta)/(2*delta);
788   }
789   else {
790     if ( r > g )
791       h = 240 + (120*(r-g)+delta)/(2*delta);
792     else
793       h = 180 + (120*(r-g+delta)+delta)/(2*delta);
794   }
795
796   return h + 1;
797 }
798
799 //-----------------------------------------------------------------------------
800 struct TVertex
801 {
802   GLfloat r, g, b, hue;
803   GLfloat vx, vy, vz;
804 };
805
806
807 //-----------------------------------------------------------------------------
808 struct TColorFunctorBase
809 {
810   virtual
811   void
812   get( TVertex& theVertex, vtkIdType thePointId, vtkIdType theCellId ) = 0;
813 };
814
815
816 //-----------------------------------------------------------------------------
817 struct TPropertyColor : TColorFunctorBase
818 {
819   vtkFloatingPointType myColor[3];
820   vtkFloatingPointType myHue;
821
822   TPropertyColor( vtkProperty *theProperty )
823   {
824     theProperty->GetColor( myColor );
825     int aRed = int( myColor[0] * 255 );
826     int aGreen = int( myColor[1] * 255 );
827     int aBlue = int( myColor[2] * 255 );
828
829     myHue = ComputeHue( aRed, aGreen, aBlue );
830   }
831
832   virtual
833   void
834   get( TVertex& theVertex, vtkIdType thePointId, vtkIdType theCellId )
835   {
836     theVertex.r = myColor[0];
837     theVertex.g = myColor[1];
838     theVertex.b = myColor[2];
839
840     theVertex.hue = myHue;
841   }
842 };
843
844
845 //-----------------------------------------------------------------------------
846 struct TColors2Color : TColorFunctorBase
847 {
848   vtkUnsignedCharArray *myColors;
849
850   TColors2Color( vtkUnsignedCharArray *theColors ):
851     myColors( theColors )
852   {}
853
854   virtual
855   void
856   get( TVertex& theVertex, vtkIdType thePointId, vtkIdType theCellId )
857   {
858     vtkIdType aTupleId = GetTupleId( thePointId, theCellId );
859     unsigned char *aColor = myColors->GetPointer( aTupleId << 2 );
860
861     theVertex.r = int( aColor[0] ) / 255.0;
862     theVertex.g = int( aColor[1] ) / 255.0;
863     theVertex.b = int( aColor[2] ) / 255.0;
864
865     theVertex.hue = ComputeHue( aColor[0], aColor[1], aColor[2] );
866   }  
867
868   virtual
869   vtkIdType
870   GetTupleId( vtkIdType thePointId, vtkIdType theCellId ) = 0;
871 };
872
873
874 //-----------------------------------------------------------------------------
875 struct TPointColors2Color : TColors2Color
876 {
877   TPointColors2Color( vtkUnsignedCharArray *theColors ):
878     TColors2Color( theColors )
879   {}
880
881   virtual
882   vtkIdType
883   GetTupleId( vtkIdType thePointId, vtkIdType theCellId )
884   {
885     return thePointId;
886   }
887 };
888
889
890 //-----------------------------------------------------------------------------
891 struct TCellColors2Color : TColors2Color
892 {
893   TCellColors2Color( vtkUnsignedCharArray *theColors ):
894     TColors2Color( theColors )
895   {}
896
897   virtual
898   vtkIdType
899   GetTupleId( vtkIdType thePointId, vtkIdType theCellId )
900   {
901     return theCellId;
902   }
903 };
904
905
906 //-----------------------------------------------------------------------------
907 template < class TCoordinates >
908 void DrawPoints( TCoordinates *theStartPoints,
909                  vtkCellArray *theCells,
910                  TColorFunctorBase* theColorFunctor,
911                  TVertex* theVertexArr,
912                  vtkIdType &theCellId,
913                  vtkIdType &theVertexId )
914 {
915   vtkIdType *ptIds = theCells->GetPointer();
916   vtkIdType *endPtIds = ptIds + theCells->GetNumberOfConnectivityEntries();
917
918   while ( ptIds < endPtIds ) {
919     vtkIdType nPts = *ptIds;
920     ++ptIds;
921
922     while ( nPts > 0 ) {
923       TVertex& aVertex = theVertexArr[ theVertexId ];
924       vtkIdType aPointId = *ptIds;
925
926       TCoordinates *anOffsetPoints = theStartPoints + 3 * aPointId;
927       aVertex.vx = anOffsetPoints[0];
928       aVertex.vy = anOffsetPoints[1];
929       aVertex.vz = anOffsetPoints[2];
930
931       theColorFunctor->get( aVertex, aPointId, theCellId );
932
933       ++theVertexId;
934       ++ptIds; 
935       --nPts; 
936     }
937
938     ++theCellId;
939   }
940 }
941
942
943 //-----------------------------------------------------------------------------
944 template < class TCoordinates >
945 void DrawCellsPoints( vtkPolyData *theInput,
946                       vtkPoints* thePoints,
947                       TColorFunctorBase* theColorFunctor,
948                       TVertex* theVertexArr )
949 {
950   vtkIdType aCellId = 0, aVertexId = 0;
951
952   TCoordinates *aStartPoints = (TCoordinates *) thePoints->GetVoidPointer(0);
953
954   if ( vtkCellArray* aCellArray = theInput->GetVerts() )
955     DrawPoints( aStartPoints, aCellArray, theColorFunctor, theVertexArr, aCellId, aVertexId );
956   
957   if ( vtkCellArray* aCellArray = theInput->GetLines() )
958     DrawPoints( aStartPoints, aCellArray, theColorFunctor, theVertexArr, aCellId, aVertexId );
959   
960   if ( vtkCellArray* aCellArray = theInput->GetPolys() )
961     DrawPoints( aStartPoints, aCellArray, theColorFunctor, theVertexArr, aCellId, aVertexId );
962   
963   if ( vtkCellArray* aCellArray = theInput->GetStrips() )
964     DrawPoints( aStartPoints, aCellArray, theColorFunctor, theVertexArr, aCellId, aVertexId ); 
965 }
966
967
968 //-----------------------------------------------------------------------------
969 int VISU_OpenGLPointSpriteMapper::Draw(vtkRenderer *theRenderer, vtkActor *theActor)
970 {
971
972   if( this->PrimitiveType == VISU_OpenGLPointSpriteMapper::GeomSphere )
973     return MAPPER_SUPERCLASS::Draw( theRenderer, theActor );
974
975   vtkUnsignedCharArray *colors = NULL;
976   vtkPolyData          *input  = this->GetInput();
977   vtkPoints            *points;
978   int noAbort = 1;
979   int cellScalars = 0;
980
981   // if the primitives are invisable then get out of here
982   if (this->PointSpriteOpacity <= 0.0)
983   {
984     return noAbort;
985   }
986
987   // and draw the display list
988   points = input->GetPoints();
989
990   // are they cell or point scalars
991   if ( this->Colors )
992   {
993     colors = this->Colors;
994     if ( (this->ScalarMode == VTK_SCALAR_MODE_USE_CELL_DATA ||
995           this->ScalarMode == VTK_SCALAR_MODE_USE_CELL_FIELD_DATA ||
996           !input->GetPointData()->GetScalars() )
997          && this->ScalarMode != VTK_SCALAR_MODE_USE_POINT_FIELD_DATA)
998     {
999       cellScalars = 1;
1000     }
1001   }
1002
1003   {
1004     vtkIdType aTotalConnectivitySize = 0;
1005
1006     if ( vtkCellArray* aCellArray = input->GetVerts() )
1007       aTotalConnectivitySize += aCellArray->GetNumberOfConnectivityEntries() - aCellArray->GetNumberOfCells();
1008
1009     if ( vtkCellArray* aCellArray = input->GetLines() )
1010       aTotalConnectivitySize += aCellArray->GetNumberOfConnectivityEntries() - aCellArray->GetNumberOfCells();
1011
1012     if ( vtkCellArray* aCellArray = input->GetPolys() )
1013       aTotalConnectivitySize += aCellArray->GetNumberOfConnectivityEntries() - aCellArray->GetNumberOfCells();
1014
1015     if ( vtkCellArray* aCellArray = input->GetStrips() )
1016       aTotalConnectivitySize += aCellArray->GetNumberOfConnectivityEntries() - aCellArray->GetNumberOfCells();
1017
1018     if ( aTotalConnectivitySize > 0 ) {
1019       TVertex* aVertexArr = new TVertex[ aTotalConnectivitySize ];
1020
1021       vtkFloatingPointType aPropertyColor[3];
1022       theActor->GetProperty()->GetColor( aPropertyColor );
1023
1024       glPointSize( this->DefaultPointSize );
1025
1026       {
1027         TColorFunctorBase* aColorFunctor = NULL;
1028         if( colors && this->PointSpriteMode != 1 ) {
1029           if ( cellScalars )
1030             aColorFunctor = new TCellColors2Color( colors );
1031           else
1032             aColorFunctor = new TPointColors2Color( colors );
1033         } else {
1034           aColorFunctor = new TPropertyColor( theActor->GetProperty() );
1035         }
1036         if ( points->GetDataType() == VTK_FLOAT )
1037           ::DrawCellsPoints< float >( input, points, aColorFunctor, aVertexArr );
1038         else
1039           ::DrawCellsPoints< double >( input, points, aColorFunctor, aVertexArr );
1040
1041         delete aColorFunctor;
1042       }
1043
1044       if( this->ExtensionsInitialized ) {
1045         GLuint aBufferObjectID = 0;
1046         vglGenBuffersARB( 1, &aBufferObjectID );
1047         vglBindBufferARB( GL_ARRAY_BUFFER_ARB, aBufferObjectID );
1048         
1049         int anArrayObjectSize = sizeof( TVertex ) * aTotalConnectivitySize;
1050         vglBufferDataARB( GL_ARRAY_BUFFER_ARB, anArrayObjectSize, aVertexArr, GL_STATIC_DRAW_ARB );
1051         
1052         delete [] aVertexArr;
1053         
1054         vglBindBufferARB( GL_ARRAY_BUFFER_ARB, 0 );
1055         vglBindBufferARB( GL_ARRAY_BUFFER_ARB, aBufferObjectID );
1056         
1057         glColorPointer( 4, GL_FLOAT, sizeof(TVertex), (void*)0 );
1058         glVertexPointer( 3, GL_FLOAT, sizeof(TVertex), (void*)(4*sizeof(GLfloat)) );
1059         
1060         glEnableClientState( GL_VERTEX_ARRAY );
1061         glEnableClientState( GL_COLOR_ARRAY );
1062         
1063         glDrawArrays( GL_POINTS, 0, aTotalConnectivitySize );
1064         
1065         glDisableClientState( GL_COLOR_ARRAY );
1066         glDisableClientState( GL_VERTEX_ARRAY );
1067         
1068         vglDeleteBuffersARB( 1, &aBufferObjectID );
1069       } else { // there are not extensions
1070         glColorPointer( 4, GL_FLOAT, sizeof(TVertex), aVertexArr );
1071         glVertexPointer( 3, GL_FLOAT, sizeof(TVertex), 
1072                          (void*)((GLfloat*)((void*)(aVertexArr)) + 4));
1073
1074         glEnableClientState( GL_VERTEX_ARRAY );
1075         glEnableClientState( GL_COLOR_ARRAY );
1076         
1077         glDrawArrays( GL_POINTS, 0, aTotalConnectivitySize );
1078         
1079         glDisableClientState( GL_COLOR_ARRAY );
1080         glDisableClientState( GL_VERTEX_ARRAY );
1081
1082         delete [] aVertexArr;
1083       }
1084     }
1085
1086     input->GetVerts()->GetNumberOfCells() + 
1087     input->GetLines()->GetNumberOfCells() + 
1088     input->GetPolys()->GetNumberOfCells() + 
1089     input->GetStrips()->GetNumberOfCells();
1090   }
1091
1092
1093   this->UpdateProgress(1.0);
1094   return noAbort;
1095 }
1096 //-----------------------------------------------------------------------------
1097 // Release the graphics resources used by this mapper.  In this case, release
1098 // the display list if any.
1099 void VISU_OpenGLPointSpriteMapper::ReleaseGraphicsResources(vtkWindow *win)
1100 {
1101   this->Superclass::ReleaseGraphicsResources(win);
1102
1103   if (this->ListId && win)
1104     {
1105     win->MakeCurrent();
1106     glDeleteLists(this->ListId,1);
1107     this->ListId = 0;
1108     }
1109   this->LastWindow = NULL;
1110 }