00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #include <config.h>
00021
00022 #ifdef HAVE_GUI_LIBS
00023
00024 #include <MeshOfHexahedra.hpp>
00025 #include <MeshOfTetrahedra.hpp>
00026 #include <Structured3DMesh.hpp>
00027
00028 #include <SurfaceMeshOfTriangles.hpp>
00029 #include <SurfaceMeshOfQuadrangles.hpp>
00030
00031 #include <RunningOptions.hpp>
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045 #include <VTKDriver.hpp>
00046
00047 #include <vtkProperty.h>
00048 #include <vtkPointData.h>
00049 #include <vtkContourFilter.h>
00050 #include <vtkCellData.h>
00051 #include <vtkLODActor.h>
00052 #include <vtkRenderer.h>
00053 #include <vtkStructuredGridOutlineFilter.h>
00054 #include <vtkPolyDataMapper.h>
00055 #include <vtkRenderWindow.h>
00056 #include <vtkXOpenGLRenderWindow.h>
00057 #include <vtkRenderWindowInteractor.h>
00058 #include <vtkTriangle.h>
00059 #include <vtkQuad.h>
00060 #include <vtkHexahedron.h>
00061 #include <vtkTetra.h>
00062 #include <vtkUnstructuredGrid.h>
00063 #include <vtkDataSetMapper.h>
00064 #include <vtkDataSetToPolyDataFilter.h>
00065 #include <vtkDoubleArray.h>
00066
00067
00068
00069
00070
00071
00072
00073
00074 vtkDoubleArray* VTKDriver::__getValues(const Mesh& m,
00075 ConstReferenceCounting<ScalarFunctionBase> pU)
00076 {
00077 vtkDoubleArray* values = vtkDoubleArray::New();
00078 values->Resize(m.numberOfVertices());
00079
00080 if (pU == 0) {
00081 for (size_t i=0; i<m.numberOfVertices(); ++i) {
00082 values->InsertNextValue(m.vertex(i).reference());
00083 }
00084 } else {
00085 #warning COULD BE OPTIMIZED A LOT
00086
00087
00088
00089
00090 const ScalarFunctionBase& u = *pU;
00091 for (size_t i=0; i<m.numberOfVertices(); ++i) {
00092 values->InsertNextValue(u(m.vertex(i)));
00093 }
00094 }
00095 return values;
00096 }
00097
00098 template <>
00099 struct VTKDriver::Traits<MeshOfTetrahedra>
00100 {
00101 typedef vtkUnstructuredGrid VTKMeshType;
00102 typedef vtkTetra VTKCellType;
00103 };
00104
00105 template <>
00106 struct VTKDriver::Traits<MeshOfHexahedra>
00107 {
00108 typedef vtkUnstructuredGrid VTKMeshType;
00109 typedef vtkHexahedron VTKCellType;
00110 };
00111
00112 template <>
00113 struct VTKDriver::Traits<Structured3DMesh>
00114 {
00115 typedef vtkUnstructuredGrid VTKMeshType;
00116 typedef vtkHexahedron VTKCellType;
00117 };
00118
00119 template <>
00120 struct VTKDriver::Traits<SurfaceMeshOfTriangles>
00121 {
00122 typedef vtkUnstructuredGrid VTKMeshType;
00123 typedef vtkTriangle VTKCellType;
00124 };
00125
00126 template <>
00127 struct VTKDriver::Traits<SurfaceMeshOfQuadrangles>
00128 {
00129 typedef vtkUnstructuredGrid VTKMeshType;
00130 typedef vtkQuad VTKCellType;
00131 };
00132
00133 template <typename MeshType>
00134 void VTKDriver::__plot(const MeshType& m,
00135 ConstReferenceCounting<ScalarFunctionBase> pU)
00136 {
00137
00138 typedef
00139 typename VTKDriver::Traits<MeshType>::VTKMeshType
00140 VTKMeshType;
00141 typedef
00142 typename VTKDriver::Traits<MeshType>::VTKCellType
00143 VTKCellType;
00144 typedef typename MeshType::CellType CellType;
00145
00146 vtkPoints* vertices = vtkPoints::New();
00147 vertices->SetNumberOfPoints(m.numberOfVertices());
00148
00149
00150 for (size_t i=0; i<m.numberOfVertices(); ++i) {
00151 const Vertex& X = m.vertex(i);
00152 vertices->InsertPoint(i,X[0],X[1],X[2]);
00153 }
00154
00155 VTKCellType* cell = VTKCellType::New();
00156 VTKMeshType* grid = VTKMeshType::New();
00157 grid->Allocate(m.numberOfCells(),1);
00158
00159
00160 for (size_t i=0; i<m.numberOfCells(); ++i) {
00161 const CellType& K = m.cell(i);
00162 for (size_t j=0; j<CellType::NumberOfVertices; ++j) {
00163 cell->GetPointIds()->SetId(j,m.vertexNumber(K(j)));
00164 }
00165 grid->InsertNextCell(cell->GetCellType(),
00166 cell->GetPointIds());
00167 }
00168
00169 grid->SetPoints(vertices);
00170
00171 vtkDoubleArray* values = this->__getValues(m, pU);
00172
00173 grid->GetPointData()->SetScalars(values);
00174
00175
00176
00177
00178
00179
00180
00181
00182 vtkDataSetMapper* mapper = vtkDataSetMapper::New();
00183 mapper->SetInput(grid);
00184 mapper->SetScalarRange(values->GetRange());
00185
00186
00187 vtkContourFilter* contour = vtkContourFilter::New();
00188 contour->SetInput(grid);
00189
00190 contour->SetNumberOfContours(1);
00191
00192
00193
00194 contour->SetValue(0,
00195 0.95*values->GetRange()[0]+0.05*values->GetRange()[1]);
00196 contour->UseScalarTreeOn();
00197
00198 vtkPolyDataMapper* mapper2 = vtkPolyDataMapper::New();
00199 mapper2->SetInput(contour->GetOutput());
00200 mapper2->SetScalarRange(values->GetRange());
00201
00202
00203
00204
00205
00206
00207 vtkLODActor *actor1 = vtkLODActor::New();
00208 actor1->SetMapper( mapper );
00209 actor1->GetProperty()->SetEdgeColor (1,1,1);
00210 actor1->GetProperty()->EdgeVisibilityOn();
00211 actor1->GetProperty()->SetRepresentationToWireframe();
00212 actor1->GetProperty()->SetOpacity(0.2);
00213 actor1->GetProperty()->SetColor(1,1,1);
00214 actor1->SetNumberOfCloudPoints(1000);
00215
00216 vtkLODActor *actor2 = vtkLODActor::New();
00217 actor2->SetMapper( mapper2 );
00218 actor2->SetNumberOfCloudPoints(1000);
00219
00220
00221
00222
00223
00224
00225 vtkRenderer *ren1= vtkRenderer::New();
00226 ren1->AddActor( actor1 );
00227 ren1->AddActor( actor2 );
00228 ren1->SetBackground( 0.1, 0.2, 0.4 );
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276 values->Delete();
00277 grid->Delete();
00278 vertices->Delete();
00279 }
00280
00281 void VTKDriver::plot(const Mesh& m,
00282 ConstReferenceCounting<ScalarFunctionBase> u)
00283 {
00284 if (not(RunningOptions::instance().useGUI())) {
00285 if (not(RunningOptions::instance().haveDisplay())) {
00286 fferr(1) << "warning: could not open DISPLAY, not processing plot()\n";
00287 } else {
00288 fferr(1) << "warning: not processing plot() when using '-nw' option\n";
00289 }
00290 return;
00291 }
00292
00293 switch (m.type()) {
00294 case Mesh::hexahedraMesh: {
00295 this->__plot(dynamic_cast<const MeshOfHexahedra&>(m), u);
00296 break;
00297 }
00298 case Mesh::cartesianHexahedraMesh: {
00299 this->__plot(dynamic_cast<const Structured3DMesh&>(m), u);
00300 break;
00301 }
00302 case Mesh::tetrahedraMesh: {
00303 this->__plot(dynamic_cast<const MeshOfTetrahedra&>(m), u);
00304 break;
00305 }
00306 case Mesh::surfaceMeshTriangles: {
00307 this->__plot(dynamic_cast<const SurfaceMeshOfTriangles&>(m), u);
00308 break;
00309 }
00310 case Mesh::surfaceMeshQuadrangles: {
00311 this->__plot(dynamic_cast<const SurfaceMeshOfQuadrangles&>(m), u);
00312 break;
00313 }
00314 default: {
00315 throw ErrorHandler(__FILE__,__LINE__,
00316 "unknown mesh type",
00317 ErrorHandler::unexpected);
00318 }
00319 }
00320 }
00321
00322 #endif // HAVE_GUI_LIBS