Changeset 3625


Ignore:
Timestamp:
Sep 19, 2006, 4:55:48 PM (18 years ago)
Author:
jack
Message:

Updatd the offline visualiser example and prepared the realtime visualiser.

Location:
anuga_core/source/anuga
Files:
2 edited
1 moved

Legend:

Unmodified
Added
Removed
  • anuga_core/source/anuga/visualiser/realtime.py

    r3549 r3625  
     1from Numeric import Float, zeros
    12from Tkinter import Button, E, W
    23from threading import Event
    34from visualiser import Visualiser
     5from vtk import vtkCellArray, vtkPoints, vtkPolyData
    46
    57class RealtimeVisualiser(Visualiser):
     
    2224
    2325    def run(self):
     26        self.tk_root.after(100, self.sync_idle.set)
    2427        Visualiser.run(self)
    25         self.tk_root.after(100, self.sync_idle.set)
     28
     29    def setup_grid(self):
     30        self.vtk_cells = vtkCellArray()
     31        # Also build vert_index - a list of the x & y values of each vertex
     32        N_vert = len(self.source.vertex_coordinates)
     33        self.vert_index = zeros((N_vert,2), Float)
     34        for n in range(N_vert):
     35            self.vtk_cells.InsertNextCell(3)
     36            for v in range(3):
     37                self.vert_index[self.source.triangles[n][v]] = self.source.vertex_coordinates[n][i*2:i*2+2]
     38                self.vtk_cells.InsertCellPoint(self.source.triangles[n][v])
     39
     40    def update_height_quantity(self, quantityName, dynamic=True):
     41        N_vert = len(self.source.vertex_coordinates)
     42        qty_index = zeros(N_vert, Float)
     43
     44        for n in range(len(self.source.triangles)):
     45            for v in range(3):
     46                qty_index[self.source.triangles[n][v]] = self.source.quantities[quantityName].vertex_values[n][v]
     47
     48        points = vtkPoints()
     49        for v in range(N_vert):
     50            points.InsertNextPoint(self.vert_index[v][0],
     51                                   self.vert_index[v][1],
     52                                   qty_index[v] * self.height_zScales[quantityName]
     53                                   + self.height_offset[quantityName])
     54
     55        polydata = self.vtk_polyData[quantityName] = vtkPolyData()
     56        polydata.SetPoints(points)
     57        polydata.SetPolys(self.vtk_cells)
     58       
     59    def build_quantity_dict(self):
     60        N_vert = len(self.source.vertex_coordinates)
     61        quantities = {}
     62        for q in self.source.quantities.keys():
     63            quantities[q] = zeros(N_vert, Float)
     64            for n in range(len(self.source.triangles)):
     65                for v in range(3):
     66                    quantities[q][self.source.triangles[n][v]] = self.source.quantities[q].vertex_values[n][v]
     67        return quantities
    2668
    2769    def setup_gui(self):
     
    4385        self.sync_idle.set()
    4486        self.sync_unpaused.set()
     87
  • anuga_core/source/anuga/visualiser/visualiser.py

    r3549 r3625  
    5151        """Create the vtkCellArray instance that represents the
    5252        triangles. Subclasses are expected to override this function
    53         to read from their source as appropriate.
     53        to read from their source as appropriate. The vtkCellArray should
     54        be stored to self.vtk_cells.
    5455        """
    5556        pass
     
    6970    def update_height_quantity(self, quantityName, dynamic=True):
    7071        """Create a vtkPolyData object and store it in
    71         self.vtk_polyData[q]. Subclasses are expected to override this
     72        self.vtk_polyData[quantityName]. Subclasses are expected to override this
    7273        function.
    7374        """
     
    115116
    116117    def build_quantity_dict(self):
    117         """Build a dictionary mapping quantity name->Numeric array of vertex
     118        """Build and return a dictionary mapping quantity name->Numeric array of vertex
    118119        values for that quantity. Subclasses are expected to override
    119120        this function."""
Note: See TracChangeset for help on using the changeset viewer.