source: inundation/ga/storm_surge/pyvolution/least_squares.py @ 390

Last change on this file since 390 was 389, checked in by duncan, 21 years ago

add vert atts to pmesh2domain/pyvolution

File size: 13.5 KB
Line 
1"""Least squares smooting and interpolation.
2
3   Implements a penalised least-squares fit and associated interpolations.
4
5   The panalty term (or smoothing term) is controlled by the smoothing
6   parameter alpha.
7   With a value of alpha=0, the fit function will attempt
8   to interpolate as closely as possible in the least-squares sense.
9   With values alpha > 0, a certain amount of smoothing will be applied.
10   A positive alpha is essential in cases where there are too few
11   data points.
12   A negative alpha is not allowed.
13   A typical value of alpha is 1.0e-6
14         
15
16   Ole Nielsen, Stephen Roberts, Duncan Gray, Christopher Zoppou
17   Geoscience Australia, 2004.   
18"""
19
20#FIXME: Current implementation uses full matrices and a general solver.
21#Later on we may consider using sparse techniques
22
23from general_mesh import General_Mesh
24from Numeric import zeros, array, Float, Int, dot, transpose
25from LinearAlgebra import solve_linear_equations
26
27try:
28    from util import gradient
29
30except ImportError, e: 
31    #FIXME reduce the dependency of modules in pyvolution
32    # Have util in a dir, working like load_mesh, and get rid of this
33    def gradient(x0, y0, x1, y1, x2, y2, q0, q1, q2):
34        """
35        """
36   
37        det = (y2-y0)*(x1-x0) - (y1-y0)*(x2-x0)           
38        a = (y2-y0)*(q1-q0) - (y1-y0)*(q2-q0)
39        a /= det
40
41        b = (x1-x0)*(q2-q0) - (x2-x0)*(q1-q0)
42        b /= det           
43
44        return a, b
45
46ALPHA_INITIAL = 0.001
47
48def fit_to_mesh_file(mesh_file, point_file, mesh_output_file, alpha= ALPHA_INITIAL):
49    """
50    Given a mesh file (tsh) and a point attribute file (xya), fit
51    point attributes to the mesh and write a mesh file with the
52    results.
53    """
54    from load_mesh.loadASCII import mesh_file_to_mesh_dictionary, \
55                 load_xya_file, export_trianglulation_file
56    # load in the .tsh file
57    mesh_dict = mesh_file_to_mesh_dictionary(mesh_file)
58    vertex_coordinates = mesh_dict['generatedpointlist']
59    triangles = mesh_dict['generatedtrianglelist']
60       
61    # load in the .xya file
62    point_dict = load_xya_file(point_file)
63    point_coordinates = point_dict['pointlist']
64    point_attributes = point_dict['pointattributelist']
65    title_string = point_dict['title']
66    title_list = title_string.split(',') #iffy! Hard coding title delimiter
67    #print "title_list", title_list   
68    for i in range(len(title_list)):
69        title_list[i] = title_list[i].strip() 
70    #print "title_list stripped", title_list   
71    f = fit_to_mesh(vertex_coordinates,
72                    triangles,
73                    point_coordinates,
74                    point_attributes,
75                    alpha = alpha)
76    # convert array to list of lists
77    mesh_dict['generatedpointattributelist'] = f.tolist()
78    mesh_dict['generatedpointattributetitlelist'] = title_list
79    export_trianglulation_file(mesh_output_file, mesh_dict)
80       
81
82def fit_to_mesh(vertex_coordinates,
83                              triangles,
84                              point_coordinates,
85                              point_attributes,
86                              alpha = ALPHA_INITIAL):
87    """
88    Fit a smooth surface to a trianglulation,
89    given data points with attributes.
90
91         
92        Inputs:
93       
94          vertex_coordinates: List of coordinate pairs [xi, eta] of points
95          constituting mesh (or a an m x 2 Numeric array)
96       
97          triangles: List of 3-tuples (or a Numeric array) of
98          integers representing indices of all vertices in the mesh.
99
100          point_coordinates: List of coordinate pairs [x, y] of data points
101          (or an nx2 Numeric array)
102
103          alpha: Smoothing parameter.
104
105          point_attributes: Vector or array of data at the point_coordinates.
106    """
107    interp = Interpolation(vertex_coordinates,
108                           triangles,
109                           point_coordinates,
110                           alpha = alpha)
111   
112    vertex_attributes = interp.fit(point_attributes)
113    return vertex_attributes
114
115
116class Interpolation:
117
118    def __init__(self,
119                 vertex_coordinates,
120                 triangles,
121                 point_coordinates = None,
122                 alpha = ALPHA_INITIAL):
123
124       
125        """ Build interpolation matrix mapping from
126        function values at vertices to function values at data points
127
128        Inputs:
129       
130          vertex_coordinates: List of coordinate pairs [xi, eta] of points
131          constituting mesh (or a an m x 2 Numeric array)
132       
133          triangles: List of 3-tuples (or a Numeric array) of
134          integers representing indices of all vertices in the mesh.
135
136          point_coordinates: List of coordinate pairs [x, y] of data points
137          (or an nx2 Numeric array)
138
139          alpha: Smoothing parameter
140         
141        """
142
143
144        #Convert input to Numeric arrays
145        vertex_coordinates = array(vertex_coordinates).astype(Float)
146        triangles = array(triangles).astype(Int)               
147       
148        #Build underlying mesh
149        self.mesh = General_Mesh(vertex_coordinates, triangles)
150
151        #Smoothing parameter
152        self.alpha = alpha
153
154        #Build coefficient matrices
155        self.build_coefficient_matrix_B(point_coordinates)   
156
157
158       
159    def build_coefficient_matrix_B(self, point_coordinates=None):
160        """Build final coefficient matrix
161        """
162
163        self.build_smoothing_matrix_D()
164       
165        if point_coordinates:
166            self.build_interpolation_matrix_A(point_coordinates)
167            AtA = dot(self.At, self.A)
168            self.B = AtA + self.alpha*self.D
169
170       
171    def build_interpolation_matrix_A(self, point_coordinates):
172        """Build n x m interpolation matrix, where
173        n is the number of data points and
174        m is the number of basis functions phi_k (one per vertex)
175        """ 
176
177        #Convert input to Numeric arrays
178        point_coordinates = array(point_coordinates).astype(Float)
179       
180        #Build n x m interpolation matrix       
181        m = self.mesh.coordinates.shape[0] #Nbr of basis functions (1/vertex)
182        n = point_coordinates.shape[0]     #Nbr of data points         
183
184        self.A = zeros((n,m), Float)
185
186        #Compute matrix elements
187        for i in range(n):
188            #For each data_coordinate point
189
190            x = point_coordinates[i]
191            for k in range(len(self.mesh)):
192                #For each triangle (brute force)
193                #FIXME: Real algorithm should only visit relevant triangles
194
195                #Get the three vertex_points
196                xi0 = self.mesh.get_vertex_coordinate(k, 0)
197                xi1 = self.mesh.get_vertex_coordinate(k, 1)
198                xi2 = self.mesh.get_vertex_coordinate(k, 2)                 
199
200                #Get the three normals
201                n0 = self.mesh.get_normal(k, 0)
202                n1 = self.mesh.get_normal(k, 1)
203                n2 = self.mesh.get_normal(k, 2)               
204
205                #Compute interpolation
206                sigma2 = dot((x-xi0), n2)/dot((xi2-xi0), n2)
207                sigma0 = dot((x-xi1), n0)/dot((xi0-xi1), n0)
208                sigma1 = dot((x-xi2), n1)/dot((xi1-xi2), n1)
209
210                #FIXME: Maybe move out to test or something
211                epsilon = 1.0e-6
212                assert abs(sigma0 + sigma1 + sigma2 - 1.0) < epsilon
213
214                #Check that this triangle contains data point
215                if sigma0 >= 0 and sigma1 >= 0 and sigma2 >= 0:
216
217                    #Assign values to matrix A
218                    j = self.mesh.triangles[k,0] #Global vertex id
219                    self.A[i, j] = sigma0
220
221                    j = self.mesh.triangles[k,1] #Global vertex id
222                    self.A[i, j] = sigma1
223
224                    j = self.mesh.triangles[k,2] #Global vertex id
225                    self.A[i, j] = sigma2
226
227   
228        #Precompute
229        self.At = transpose(self.A)
230
231
232        #FIXME: Remember to re-introduce the 1/n factor in the
233        #interpolation term
234       
235    def build_smoothing_matrix_D(self):
236        """Build m x m smoothing matrix, where
237        m is the number of basis functions phi_k (one per vertex)
238
239        The smoothing matrix is defined as
240
241        D = D1 + D2
242
243        where
244
245        [D1]_{k,l} = \int_\Omega
246           \frac{\partial \phi_k}{\partial x}
247           \frac{\partial \phi_l}{\partial x}\,
248           dx dy
249
250        [D2]_{k,l} = \int_\Omega
251           \frac{\partial \phi_k}{\partial y}
252           \frac{\partial \phi_l}{\partial y}\,
253           dx dy
254
255
256        The derivatives \frac{\partial \phi_k}{\partial x},
257        \frac{\partial \phi_k}{\partial x} for a particular triangle
258        are obtained by computing the gradient a_k, b_k for basis function k
259        """
260
261        #FIXME: algorithm might be optimised by computing local 9x9
262        #"element stiffness matrices:
263
264       
265        m = self.mesh.coordinates.shape[0] #Nbr of basis functions (1/vertex)
266
267        self.D = zeros((m,m), Float)
268
269        #For each triangle compute contributions to D = D1+D2       
270        for i in range(len(self.mesh)):
271
272            #Get area
273            area = self.mesh.areas[i]
274
275            #Get global vertex indices
276            v0 = self.mesh.triangles[i,0]
277            v1 = self.mesh.triangles[i,1]
278            v2 = self.mesh.triangles[i,2]
279           
280            #Get the three vertex_points
281            xi0 = self.mesh.get_vertex_coordinate(i, 0)
282            xi1 = self.mesh.get_vertex_coordinate(i, 1)
283            xi2 = self.mesh.get_vertex_coordinate(i, 2)                 
284
285            #Compute gradients for each vertex
286            a0, b0 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
287                              1, 0, 0)
288
289            a1, b1 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
290                              0, 1, 0)
291
292            a2, b2 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
293                              0, 0, 1)           
294
295            #Compute diagonal contributions
296            self.D[v0,v0] += (a0*a0 + b0*b0)*area
297            self.D[v1,v1] += (a1*a1 + b1*b1)*area
298            self.D[v2,v2] += (a2*a2 + b2*b2)*area           
299
300            #Compute contributions for basis functions sharing edges
301            e01 = (a0*a1 + b0*b1)*area
302            self.D[v0,v1] += e01
303            self.D[v1,v0] += e01
304
305            e12 = (a1*a2 + b1*b2)*area
306            self.D[v1,v2] += e12
307            self.D[v2,v1] += e12
308
309            e20 = (a2*a0 + b2*b0)*area
310            self.D[v2,v0] += e20
311            self.D[v0,v2] += e20             
312           
313           
314    def fit(self, z):
315        """Fit a smooth surface to given data points z.
316
317        The smooth surface is computed at each vertex in the underlying
318        mesh using the formula given in the module doc string.
319
320        Pre Condition:
321          self.A, self.At and self.B have been initialised
322         
323        Inputs:
324          z: Vector or array of data at the point_coordinates.
325          If z is an array, smoothing will be done for each column
326        """
327
328        #Convert input to Numeric arrays
329        z = array(z).astype(Float)
330       
331        #Compute right hand side based on data
332        Atz = dot(self.At, z)
333
334        #Check sanity
335        n, m = self.A.shape
336        if n<m and self.alpha == 0.0:
337            msg = 'ERROR (least_squares): Too few data points\n'
338            msg += 'There only %d data points. Need at least %d\n' %(n,m)
339            msg += 'Alternatively, increase smoothing parameter alpha' 
340            raise msg
341
342
343        #Solve and return
344        return solve_linear_equations(self.B, Atz)
345
346        #FIXME: Should we store the result here for later use? (ON)
347   
348
349    def interpolate(self, f):
350        """Compute predicted values at data points implied in self.A.
351
352        The mesh values representing a smooth surface are
353        assumed to be specified in f. This argument could,
354        for example have been obtained from the method self.fit()
355       
356        Pre Condition:
357          self.A has been initialised
358       
359        Inputs:
360          f: Vector or array of data at the mesh vertices.
361          If f is an array, interpolation will be done for each column
362        """
363       
364        return dot(self.A, f)
365
366     
367    #FIXME: We will need a method 'evaluate(self):' that will interpolate
368    #a computed surface living on the mesh onto a collection of
369    #arbitrary data points
370    #
371    #Precondition: self.fit(z) has stored its result in self.f.
372    #
373    #Input: data_points
374    #Algorithm:
375    #  1 Build a new temporary A matrix based on mesh and new data points
376    #  2 Apply it to self.f (return A*self.f)
377    #
378    # ON
379
380
381
382#-------------------------------------------------------------
383if __name__ == "__main__":
384    """
385    Load in a mesh and data points with attributes.
386    Fit the attributes to the mesh.
387    Save a new mesh file.
388    """
389    import os, sys
390    usage = "usage: %s mesh_input.tsh point.xya  mesh_output.tsh alpha" %         os.path.basename(sys.argv[0])
391
392    if len(sys.argv) < 4:
393        print usage
394    else:
395        mesh_file = sys.argv[1]
396        point_file = sys.argv[2]
397        mesh_output_file = sys.argv[3]
398        if len(sys.argv) > 4:
399            alpha = sys.argv[4]
400        else:
401            alpha = ALPHA_INITIAL
402        fit_to_mesh_file(mesh_file, point_file, mesh_output_file, alpha)
403       
Note: See TracBrowser for help on using the repository browser.