yade-dev team mailing list archive
-
yade-dev team
-
Mailing list archive
-
Message #02028
[Branch ~yade-dev/yade/trunk] Rev 1758: Added a c++ VTKRecorder (for spheres and facets) and new feature 'vtk' needed for it.
------------------------------------------------------------
revno: 1758
committer: Sergei D. <sega@laptop>
branch nick: trunk
timestamp: Tue 2009-09-15 23:32:10 +0400
message:
Added a c++ VTKRecorder (for spheres and facets) and new feature 'vtk' needed for it.
added:
pkg/dem/Engine/StandAloneEngine/VTKRecorder.cpp
pkg/dem/Engine/StandAloneEngine/VTKRecorder.hpp
modified:
SConstruct
py/export.py
--
lp:yade
https://code.launchpad.net/~yade-dev/yade/trunk
Your team Yade developers is subscribed to branch lp:yade.
To unsubscribe from this branch go to https://code.launchpad.net/~yade-dev/yade/trunk/+edit-subscription.
=== modified file 'SConstruct'
--- SConstruct 2009-08-28 13:23:39 +0000
+++ SConstruct 2009-09-15 19:32:10 +0000
@@ -136,7 +136,7 @@
ListVariable('exclude','Yade components that will not be built','none',names=['qt3','gui','extra','common','dem','fem','lattice','mass-spring','realtime-rigidbody','snow']),
EnumVariable('PGO','Whether to "gen"erate or "use" Profile-Guided Optimization','',['','gen','use'],{'no':'','0':'','false':''},1),
# OK, dummy prevents bug in scons: if one selects all, it says all in scons.config, but without quotes, which generates error.
- ListVariable('features','Optional features that are turned on','python,log4cxx,opengl,gts,openmp',names=['opengl','python','log4cxx','cgal','openmp','gts']),
+ ListVariable('features','Optional features that are turned on','python,log4cxx,opengl,gts,openmp',names=['opengl','python','log4cxx','cgal','openmp','gts','vtk']),
('jobs','Number of jobs to run at the same time (same as -j, but saved)',4,None,int),
('extraModules', 'Extra directories with their own SConscript files (must be in-tree) (whitespace separated)',None,None,Split),
('buildPrefix','Where to create build-[version][variant] directory for intermediary files','..'),
@@ -327,6 +327,9 @@
print "\nERROR: Unable to compile with optional feature `%s'.\n\nIf you are sure, remove it from features (scons features=featureOne,featureTwo for example) and build again."%featureName
Exit(1)
# check "optional" libs
+ if 'vtk' in env['features']:
+ ok=conf.CheckLibWithHeader('vtkHybrid','vtkInstantiator.h','c++','vtkInstantiator::New();',autoadd=1)
+ if not ok: featureNotOK('vtk')
if 'opengl' in env['features']:
ok=conf.CheckLibWithHeader('glut','GL/glut.h','c++','glutGetModifiers();',autoadd=1)
# TODO ok=True for darwin platform where openGL (and glut) is native
=== added file 'pkg/dem/Engine/StandAloneEngine/VTKRecorder.cpp'
--- pkg/dem/Engine/StandAloneEngine/VTKRecorder.cpp 1970-01-01 00:00:00 +0000
+++ pkg/dem/Engine/StandAloneEngine/VTKRecorder.cpp 2009-09-15 19:32:10 +0000
@@ -0,0 +1,121 @@
+#include"VTKRecorder.hpp"
+#include<vtkCellArray.h>
+#include<vtkPoints.h>
+#include<vtkPointData.h>
+#include<vtkSmartPointer.h>
+#include<vtkFloatArray.h>
+#include<vtkUnstructuredGrid.h>
+#include<vtkXMLUnstructuredGridWriter.h>
+#include<vtkXMLMultiBlockDataWriter.h>
+#include<vtkMultiBlockDataSet.h>
+#include<vtkTriangle.h>
+#include<yade/core/MetaBody.hpp>
+#include<yade/pkg-common/Sphere.hpp>
+#include<yade/pkg-common/Facet.hpp>
+
+YADE_PLUGIN((VTKRecorder));
+YADE_REQUIRE_FEATURE(VTK)
+CREATE_LOGGER(VTKRecorder);
+
+VTKRecorder::VTKRecorder()
+{
+ /* we always want to save the first state as well */
+ initRun=true;
+}
+
+VTKRecorder::~VTKRecorder()
+{
+
+}
+
+void VTKRecorder::init(MetaBody* rootBody)
+{
+}
+
+void VTKRecorder::action(MetaBody* rootBody)
+{
+ vector<bool> recActive(REC_SENTINEL,false);
+ FOREACH(string& rec, recorders){
+ if(rec=="spheres") recActive[REC_SPHERES]=true;
+ else if(rec=="facets") recActive[REC_FACETS]=true;
+ else LOG_ERROR("Unknown recorder named `"<<rec<<"' (supported are: spheres, facets). Ignored.");
+ }
+
+ vtkSmartPointer<vtkPoints> spheresPos = vtkSmartPointer<vtkPoints>::New();
+ vtkSmartPointer<vtkCellArray> spheresCells = vtkSmartPointer<vtkCellArray>::New();
+ vtkSmartPointer<vtkFloatArray> radii = vtkSmartPointer<vtkFloatArray>::New();
+ radii->SetNumberOfComponents(1);
+ radii->SetName("Radii");
+
+ vtkSmartPointer<vtkPoints> facetsPos = vtkSmartPointer<vtkPoints>::New();
+ vtkSmartPointer<vtkCellArray> facetsCells = vtkSmartPointer<vtkCellArray>::New();
+
+ FOREACH(const shared_ptr<Body>& b, *rootBody->bodies){
+ if (recActive[REC_SPHERES])
+ {
+ const Sphere* sphere = dynamic_cast<Sphere*>(b->geometricalModel.get());
+ if (sphere)
+ {
+ vtkIdType pid[1];
+ const Vector3r& pos = b->physicalParameters->se3.position;
+ pid[0] = spheresPos->InsertNextPoint(pos[0], pos[1], pos[2]);
+ spheresCells->InsertNextCell(1,pid);
+ radii->InsertNextValue(sphere->radius);
+ continue;
+ }
+ }
+ if (recActive[REC_FACETS])
+ {
+ const Facet* facet = dynamic_cast<Facet*>(b->geometricalModel.get());
+ if (facet)
+ {
+ const Se3r& O = b->physicalParameters->se3;
+ const vector<Vector3r>& localPos = facet->vertices;
+ Matrix3r facetAxisT; O.orientation.ToRotationMatrix(facetAxisT);
+ vtkSmartPointer<vtkTriangle> tri = vtkSmartPointer<vtkTriangle>::New();
+ vtkIdType nbPoints=facetsPos->GetNumberOfPoints();
+ for (int i=0;i<3;++i) {
+ Vector3r globalPos = O.position + facetAxisT * localPos[i];
+ facetsPos->InsertNextPoint(globalPos[0], globalPos[1], globalPos[2]);
+ tri->GetPointIds()->SetId(i,nbPoints+i);
+ }
+ facetsCells->InsertNextCell(tri);
+ continue;
+ }
+ }
+ }
+
+ if (recActive[REC_SPHERES])
+ {
+ vtkSmartPointer<vtkUnstructuredGrid> spheresUg = vtkSmartPointer<vtkUnstructuredGrid>::New();
+ spheresUg->SetPoints(spheresPos);
+ spheresUg->SetCells(VTK_VERTEX, spheresCells);
+ spheresUg->GetPointData()->AddArray(radii);
+ vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
+ string fn=fileName+"spheres."+lexical_cast<string>(rootBody->currentIteration)+".vtu";
+ writer->SetFileName(fn.c_str());
+ writer->SetInput(spheresUg);
+ writer->Write();
+ }
+ if (recActive[REC_FACETS])
+ {
+ vtkSmartPointer<vtkUnstructuredGrid> facetsUg = vtkSmartPointer<vtkUnstructuredGrid>::New();
+ facetsUg->SetPoints(facetsPos);
+ facetsUg->SetCells(VTK_TRIANGLE, facetsCells);
+ vtkSmartPointer<vtkXMLUnstructuredGridWriter> writer = vtkSmartPointer<vtkXMLUnstructuredGridWriter>::New();
+ string fn=fileName+"facets."+lexical_cast<string>(rootBody->currentIteration)+".vtu";
+ writer->SetFileName(fn.c_str());
+ writer->SetInput(facetsUg);
+ writer->Write();
+ }
+
+ //vtkSmartPointer<vtkMultiBlockDataSet> multiblockDataset = vtkSmartPointer<vtkMultiBlockDataSet>::New();
+ //multiblockDataset->SetBlock(0, spheresUg );
+ //multiblockDataset->SetBlock(1, facetsUg );
+ //vtkSmartPointer<vtkXMLMultiBlockDataWriter> writer = vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New();
+ //string fn=fileName+lexical_cast<string>(rootBody->currentIteration)+".vtm";
+ //writer->SetFileName(fn.c_str());
+ //writer->SetInput(multiblockDataset);
+ //writer->Write();
+}
+
=== added file 'pkg/dem/Engine/StandAloneEngine/VTKRecorder.hpp'
--- pkg/dem/Engine/StandAloneEngine/VTKRecorder.hpp 1970-01-01 00:00:00 +0000
+++ pkg/dem/Engine/StandAloneEngine/VTKRecorder.hpp 2009-09-15 19:32:10 +0000
@@ -0,0 +1,22 @@
+#pragma once
+#include<yade/pkg-common/PeriodicEngines.hpp>
+
+
+class VTKRecorder: public PeriodicEngine {
+ public:
+ enum {REC_SPHERES=0,REC_FACETS,REC_SENTINEL};
+ vector<string> recorders;
+ string fileName;
+ //bool multiBlockData;
+ VTKRecorder(); //{ /* we always want to save the first state as well */ initRun=true; };
+ ~VTKRecorder();
+ void init(MetaBody*);
+ virtual void action(MetaBody*);
+ private:
+
+ REGISTER_ATTRIBUTES(PeriodicEngine,(recorders)(fileName));
+ REGISTER_CLASS_AND_BASE(VTKRecorder,PeriodicEngine);
+ DECLARE_LOGGER;
+};
+REGISTER_SERIALIZABLE(VTKRecorder);
+
=== modified file 'py/export.py'
--- py/export.py 2009-09-05 10:11:06 +0000
+++ py/export.py 2009-09-15 19:32:10 +0000
@@ -111,7 +111,7 @@
piece.appendChild(cell_data)
# Write to file and exit
- outFile = open(self.baseName+'-%04d'%self.snapCount+'.vtu', 'w')
+ outFile = open(self.baseName+'%04d'%self.snapCount+'.vtu', 'w')
# xml.dom.ext.PrettyPrint(doc, file)
doc.writexml(outFile, newl='\n')
outFile.close()