source: anuga_core/source/anuga/visualiser/realtime.py @ 3958

Last change on this file since 3958 was 3958, checked in by jack, 17 years ago

Updated the visualisers to work with the new general_mesh structure.

File size: 5.7 KB
Line 
1from Numeric import Float, zeros
2from Tkinter import Button, E, Tk, W
3from threading import Event
4from visualiser import Visualiser
5from vtk import vtkCellArray, vtkPoints, vtkPolyData
6
7class RealtimeVisualiser(Visualiser):
8    """A VTK-powered realtime visualiser which runs in its own thread.
9    In addition to the functions provided by the standard visualiser,
10    the following additional functions are provided:
11
12    update() - Sync the visualiser to the current state of the model.
13    Should be called inside the evolve loop.
14
15    evolveFinished() - Clean up synchronisation constructs that tie the
16    visualiser to the evolve loop. Call this after the evolve loop finishes
17    to ensure a clean shutdown.
18    """
19    def __init__(self, source):
20        """The source parameter is assumed to be a Domain.
21        """
22        Visualiser.__init__(self, source)
23
24        self.running = True
25
26        self.xmin = None
27        self.xmax = None
28        self.ymin = None
29        self.ymax = None
30        self.zmin = None
31        self.zmax = None
32
33        # Synchronisation Constructs
34        self.sync_idle = Event()
35        self.sync_idle.clear()
36        self.sync_unpaused = Event()
37        self.sync_unpaused.set()
38        self.sync_redrawReady = Event()
39        self.sync_redrawReady.clear()
40
41    def run(self):
42        self.alter_tkroot(Tk.after, (100, self.sync_idle.set))
43        Visualiser.run(self)
44
45    def setup_grid(self):
46        self.vtk_cells = vtkCellArray()
47        triangles = self.source.triangles
48        N_tri = self.source.number_of_triangles
49        verticies = self.source.get_vertex_coordinates()
50        N_vert = len(verticies)
51        # Also build vert_index - a list of the x & y values of each vertex
52        self.vert_index = zeros((N_vert,2), Float)
53        for n in range(N_tri):
54            self.vtk_cells.InsertNextCell(3)
55            for v in range(3):
56                self.vert_index[triangles[n][v]] = verticies[n * 3 + v]
57                self.vtk_cells.InsertCellPoint(triangles[n][v])
58
59    def update_height_quantity(self, quantityName, dynamic=True):
60        N_vert = len(self.source.vertex_coordinates)
61        qty_index = zeros(N_vert, Float)
62
63        for n in range(len(self.source.triangles)):
64            for v in range(3):
65                qty_index[self.source.triangles[n][v]] = self.source.quantities[quantityName].vertex_values[n][v]
66
67        points = vtkPoints()
68        for v in range(N_vert):
69            points.InsertNextPoint(self.vert_index[v][0],
70                                   self.vert_index[v][1],
71                                   qty_index[v] * self.height_zScales[quantityName]
72                                   + self.height_offset[quantityName])
73            if self.xmin == None or self.xmin > self.vert_index[v][0]:
74                self.xmin = self.vert_index[v][0]
75            if self.xmax == None or self.xmax < self.vert_index[v][0]:
76                self.xmax = self.vert_index[v][0]
77            if self.ymin == None or self.ymin > self.vert_index[v][1]:
78                self.ymin = self.vert_index[v][1]
79            if self.ymax == None or self.ymax < self.vert_index[v][1]:
80                self.ymax = self.vert_index[v][1]
81            if self.zmin == None or self.zmin > qty_index[v] * self.height_zScales[quantityName] + self.height_offset[quantityName]:
82                self.zmin = qty_index[v] * self.height_zScales[quantityName] + self.height_offset[quantityName]
83            if self.zmax == None or self.zmax < qty_index[v] * self.height_zScales[quantityName] + self.height_offset[quantityName]:
84                self.zmax = qty_index[v] * self.height_zScales[quantityName] + self.height_offset[quantityName]
85
86        polydata = self.vtk_polyData[quantityName] = vtkPolyData()
87        polydata.SetPoints(points)
88        polydata.SetPolys(self.vtk_cells)
89
90    def get_3d_bounds(self):
91        return [self.xmin, self.xmax, self.ymin, self.ymax, self.zmin, self.zmax]
92       
93    def build_quantity_dict(self):
94        N_vert = len(self.source.vertex_coordinates)
95        quantities = {}
96        for q in self.source.quantities.keys():
97            quantities[q] = zeros(N_vert, Float)
98            for n in range(len(self.source.triangles)):
99                for v in range(3):
100                    quantities[q][self.source.triangles[n][v]] = self.source.quantities[q].vertex_values[n][v]
101        return quantities
102
103    def setup_gui(self):
104        Visualiser.setup_gui(self)
105        self.tk_pauseResume = Button(self.tk_controlFrame, text="Pause", command=self.pauseResume)
106        self.tk_pauseResume.grid(row=1, column=0, sticky=E+W)
107
108    def pauseResume(self):
109        if self.sync_unpaused.isSet():
110            self.sync_unpaused.clear()
111            self.tk_pauseResume.config(text="Resume")
112        else:
113            self.sync_unpaused.set()
114            self.tk_pauseResume.config(text="Pause")
115
116    def shutdown(self):
117        Visualiser.shutdown(self)
118        self.running = False
119        self.sync_idle.set()
120        self.sync_unpaused.set()
121
122    def redraw(self):
123        if self.running and self.sync_unpaused.isSet():
124            self.sync_redrawReady.wait()
125            self.sync_redrawReady.clear()
126            self.redraw_quantities()
127            self.sync_idle.set()
128        Visualiser.redraw(self)
129
130    def update(self):
131        """Sync the visualiser to the domain. Call this in the evolve loop."""
132        if self.running:
133            self.sync_redrawReady.set()
134            self.sync_idle.wait()
135            self.sync_idle.clear()
136            self.sync_unpaused.wait()
137
138    def evolveFinished(self):
139        """Stop the visualiser from waiting on signals from the evolve loop.
140        Call this just after the evolve loop to ensure a clean shutdown."""
141        self.running = False
142        self.sync_redrawReady.set()
Note: See TracBrowser for help on using the repository browser.