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

Last change on this file since 474 was 467, checked in by ole, 21 years ago
File size: 17.2 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
21#FIXME: Current implementation uses full matrices and a general solver.
22#Later on we may consider using sparse techniques
23
24import exceptions
25class ShapeError(exceptions.Exception): pass
26
27from general_mesh import General_Mesh
28from Numeric import zeros, array, Float, Int, dot, transpose
29from LinearAlgebra import solve_linear_equations
30#from scipy import sparse
31from sparse import Sparse
32from cg_solve import conjugate_gradient, VectorShapeError
33
34try:
35    from util import gradient
36
37except ImportError, e: 
38    #FIXME reduce the dependency of modules in pyvolution
39    # Have util in a dir, working like load_mesh, and get rid of this
40    def gradient(x0, y0, x1, y1, x2, y2, q0, q1, q2):
41        """
42        """
43   
44        det = (y2-y0)*(x1-x0) - (y1-y0)*(x2-x0)           
45        a = (y2-y0)*(q1-q0) - (y1-y0)*(q2-q0)
46        a /= det
47
48        b = (x1-x0)*(q2-q0) - (x2-x0)*(q1-q0)
49        b /= det           
50
51        return a, b
52
53
54DEFAULT_ALPHA = 0.001
55   
56
57def fit_to_mesh_file(mesh_file, point_file, mesh_output_file, alpha=DEFAULT_ALPHA):
58    """
59    Given a mesh file (tsh) and a point attribute file (xya), fit
60    point attributes to the mesh and write a mesh file with the
61    results.
62    """
63    from load_mesh.loadASCII import mesh_file_to_mesh_dictionary, \
64                 load_xya_file, export_trianglulation_file
65    # load in the .tsh file
66    mesh_dict = mesh_file_to_mesh_dictionary(mesh_file)
67    vertex_coordinates = mesh_dict['generatedpointlist']
68    triangles = mesh_dict['generatedtrianglelist']
69   
70    old_point_attributes = mesh_dict['generatedpointattributelist'] 
71    old_title_list = mesh_dict['generatedpointattributetitlelist']
72
73   
74    # load in the .xya file
75   
76    point_dict = load_xya_file(point_file)
77    point_coordinates = point_dict['pointlist']
78    point_attributes = point_dict['pointattributelist']
79    title_string = point_dict['title']
80    title_list = title_string.split(',') #iffy! Hard coding title delimiter 
81    for i in range(len(title_list)):
82        title_list[i] = title_list[i].strip() 
83    #print "title_list stripped", title_list   
84    f = fit_to_mesh(vertex_coordinates,
85                    triangles,
86                    point_coordinates,
87                    point_attributes,
88                    alpha = alpha)
89   
90    # convert array to list of lists
91    new_point_attributes = f.tolist()
92
93    #FIXME have this overwrite attributes with the same title - DSG
94    #Put the newer attributes last
95    if old_title_list <> []:
96        old_title_list.extend(title_list)
97        #FIXME can this be done a faster way? - DSG
98        for i in range(len(old_point_attributes)):
99            old_point_attributes[i].extend(new_point_attributes[i])
100        mesh_dict['generatedpointattributelist'] = old_point_attributes
101        mesh_dict['generatedpointattributetitlelist'] = old_title_list
102    else:
103        mesh_dict['generatedpointattributelist'] = new_point_attributes
104        mesh_dict['generatedpointattributetitlelist'] = title_list
105   
106    export_trianglulation_file(mesh_output_file, mesh_dict)
107       
108
109def fit_to_mesh(vertex_coordinates,
110                triangles,
111                point_coordinates,
112                point_attributes,
113                alpha = DEFAULT_ALPHA):
114    """
115    Fit a smooth surface to a trianglulation,
116    given data points with attributes.
117
118         
119        Inputs:
120       
121          vertex_coordinates: List of coordinate pairs [xi, eta] of points
122          constituting mesh (or a an m x 2 Numeric array)
123       
124          triangles: List of 3-tuples (or a Numeric array) of
125          integers representing indices of all vertices in the mesh.
126
127          point_coordinates: List of coordinate pairs [x, y] of data points
128          (or an nx2 Numeric array)
129
130          alpha: Smoothing parameter.
131
132          point_attributes: Vector or array of data at the point_coordinates.
133    """
134    interp = Interpolation(vertex_coordinates,
135                           triangles,
136                           point_coordinates,
137                           alpha = alpha)
138   
139    vertex_attributes = interp.fit_points(point_attributes)
140    return vertex_attributes
141
142
143class Interpolation:
144
145    def __init__(self,
146                 vertex_coordinates,
147                 triangles,
148                 point_coordinates = None,
149                 alpha = DEFAULT_ALPHA):
150
151       
152        """ Build interpolation matrix mapping from
153        function values at vertices to function values at data points
154
155        Inputs:
156       
157          vertex_coordinates: List of coordinate pairs [xi, eta] of points
158          constituting mesh (or a an m x 2 Numeric array)
159       
160          triangles: List of 3-tuples (or a Numeric array) of
161          integers representing indices of all vertices in the mesh.
162
163          point_coordinates: List of coordinate pairs [x, y] of data points
164          (or an nx2 Numeric array)
165
166          alpha: Smoothing parameter
167         
168        """
169
170
171        #Convert input to Numeric arrays
172        vertex_coordinates = array(vertex_coordinates).astype(Float)
173        triangles = array(triangles).astype(Int)               
174       
175        #Build underlying mesh
176        self.mesh = General_Mesh(vertex_coordinates, triangles)
177
178        #Smoothing parameter
179        self.alpha = alpha
180
181        #Build coefficient matrices
182        self.build_coefficient_matrix_B(point_coordinates)   
183
184
185       
186    def build_coefficient_matrix_B(self, point_coordinates=None):
187        """Build final coefficient matrix"""
188       
189
190        if self.alpha <> 0:
191            self.build_smoothing_matrix_D()
192       
193        if point_coordinates:
194
195            self.build_interpolation_matrix_A(point_coordinates)
196
197            if self.alpha <> 0:
198                self.B = self.AtA + self.alpha*self.D
199            else:
200                self.B = self.AtA
201
202       
203    def build_interpolation_matrix_A(self, point_coordinates):
204        """Build n x m interpolation matrix, where
205        n is the number of data points and
206        m is the number of basis functions phi_k (one per vertex)"""
207       
208        #Convert input to Numeric arrays
209        point_coordinates = array(point_coordinates).astype(Float)
210       
211        #Build n x m interpolation matrix       
212        m = self.mesh.coordinates.shape[0] #Nbr of basis functions (1/vertex)
213        n = point_coordinates.shape[0]     #Nbr of data points
214       
215        #self.A = zeros((n,m), Float)
216        self.A = Sparse(n,m)
217        self.AtA = Sparse(m,m)
218
219        #Compute matrix elements
220        for i in range(n):
221            #For each data_coordinate point
222
223            #print 'Doing %d/%d' %(i, n)
224
225            x = point_coordinates[i]
226            element_found = False
227            k = 0
228            while not element_found and k < len(self.mesh):
229                #For each triangle (brute force)
230                #FIXME: Real algorithm should only visit relevant triangles
231
232                #Get the three vertex_points
233                xi0 = self.mesh.get_vertex_coordinate(k, 0)
234                xi1 = self.mesh.get_vertex_coordinate(k, 1)
235                xi2 = self.mesh.get_vertex_coordinate(k, 2)                 
236
237                #Get the three normals
238                n0 = self.mesh.get_normal(k, 0)
239                n1 = self.mesh.get_normal(k, 1)
240                n2 = self.mesh.get_normal(k, 2)               
241
242                #Compute interpolation
243                sigma2 = dot((x-xi0), n2)/dot((xi2-xi0), n2)
244                sigma0 = dot((x-xi1), n0)/dot((xi0-xi1), n0)
245                sigma1 = dot((x-xi2), n1)/dot((xi1-xi2), n1)
246
247                #FIXME: Maybe move out to test or something
248                epsilon = 1.0e-6
249                assert abs(sigma0 + sigma1 + sigma2 - 1.0) < epsilon
250
251                #Check that this triangle contains data point
252                if sigma0 >= 0 and sigma1 >= 0 and sigma2 >= 0:
253                    element_found = True
254                    #Assign values to matrix A
255
256
257                    j0 = self.mesh.triangles[k,0] #Global vertex id
258                    #self.A[i, j0] = sigma0
259
260                    j1 = self.mesh.triangles[k,1] #Global vertex id
261                    #self.A[i, j1] = sigma1
262
263                    j2 = self.mesh.triangles[k,2] #Global vertex id
264                    #self.A[i, j2] = sigma2
265
266                    sigmas = {j0:sigma0, j1:sigma1, j2:sigma2}
267                    js     = [j0,j1,j2]
268
269                    for j in js:
270                        self.A[i,j] = sigmas[j]
271                        for k in js:
272                            self.AtA[j,k] += sigmas[j]*sigmas[k]
273                k = k+1
274       
275##         self.A   = (self.A).tocsc()
276##         self.AtA = (self.AtA).tocsc()
277##         self.At  =  self.A.transp()
278
279
280       
281    def get_A(self):
282        return self.A.todense() 
283
284    def get_B(self):
285        return self.B.todense()
286   
287    def get_D(self):
288        return self.D.todense()
289   
290        #FIXME: Remember to re-introduce the 1/n factor in the
291        #interpolation term
292       
293    def build_smoothing_matrix_D(self):
294        """Build m x m smoothing matrix, where
295        m is the number of basis functions phi_k (one per vertex)
296
297        The smoothing matrix is defined as
298
299        D = D1 + D2
300
301        where
302
303        [D1]_{k,l} = \int_\Omega
304           \frac{\partial \phi_k}{\partial x}
305           \frac{\partial \phi_l}{\partial x}\,
306           dx dy
307
308        [D2]_{k,l} = \int_\Omega
309           \frac{\partial \phi_k}{\partial y}
310           \frac{\partial \phi_l}{\partial y}\,
311           dx dy
312
313
314        The derivatives \frac{\partial \phi_k}{\partial x},
315        \frac{\partial \phi_k}{\partial x} for a particular triangle
316        are obtained by computing the gradient a_k, b_k for basis function k
317        """
318
319        #FIXME: algorithm might be optimised by computing local 9x9
320        #"element stiffness matrices:
321
322       
323        m = self.mesh.coordinates.shape[0] #Nbr of basis functions (1/vertex)
324
325        #self.D = zeros((m,m), Float)
326        self.D = Sparse(m,m)
327
328        #For each triangle compute contributions to D = D1+D2       
329        for i in range(len(self.mesh)):
330
331            #Get area
332            area = self.mesh.areas[i]
333
334            #Get global vertex indices
335            v0 = self.mesh.triangles[i,0]
336            v1 = self.mesh.triangles[i,1]
337            v2 = self.mesh.triangles[i,2]
338           
339            #Get the three vertex_points
340            xi0 = self.mesh.get_vertex_coordinate(i, 0)
341            xi1 = self.mesh.get_vertex_coordinate(i, 1)
342            xi2 = self.mesh.get_vertex_coordinate(i, 2)                 
343
344            #Compute gradients for each vertex
345            a0, b0 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
346                              1, 0, 0)
347
348            a1, b1 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
349                              0, 1, 0)
350
351            a2, b2 = gradient(xi0[0], xi0[1], xi1[0], xi1[1], xi2[0], xi2[1],
352                              0, 0, 1)           
353
354            #Compute diagonal contributions
355            self.D[v0,v0] += (a0*a0 + b0*b0)*area
356            self.D[v1,v1] += (a1*a1 + b1*b1)*area
357            self.D[v2,v2] += (a2*a2 + b2*b2)*area           
358
359            #Compute contributions for basis functions sharing edges
360            e01 = (a0*a1 + b0*b1)*area
361            self.D[v0,v1] += e01
362            self.D[v1,v0] += e01
363
364            e12 = (a1*a2 + b1*b2)*area
365            self.D[v1,v2] += e12
366            self.D[v2,v1] += e12
367
368            e20 = (a2*a0 + b2*b0)*area
369            self.D[v2,v0] += e20
370            self.D[v0,v2] += e20             
371
372        #self.D = (self.D).tocsc()
373           
374    def fit(self, z):
375        """Fit a smooth surface to given 1d array of data points z.
376
377        The smooth surface is computed at each vertex in the underlying
378        mesh using the formula given in the module doc string.
379
380        Pre Condition:
381          self.A, self.At and self.B have been initialised
382         
383        Inputs:
384          z: Single 1d vector or array of data at the point_coordinates.
385        """
386
387        #Convert input to Numeric arrays
388        z = array(z).astype(Float)
389
390
391        if len(z.shape) > 1 :
392            raise VectorShapeError, 'Can only deal with 1d data vector'
393       
394        #Compute right hand side based on data
395        Atz = self.A.trans_mult(z)
396
397
398        #print 'fit: Atz',Atz
399       
400        #Check sanity
401        n, m = self.A.shape
402        if n<m and self.alpha == 0.0:
403            msg = 'ERROR (least_squares): Too few data points\n'
404            msg += 'There only %d data points. Need at least %d\n' %(n,m)
405            msg += 'Alternatively, increase smoothing parameter alpha' 
406            raise msg
407
408
409        #Solve and return
410        #return solve_linear_equations(self.get_B(), Atz)
411
412        # caused errors
413        #return sparse.solve(self.B,Atz)
414
415        return conjugate_gradient(self.B, Atz, Atz)
416        #FIXME: Should we store the result here for later use? (ON)       
417   
418    def fit_points(self, z):
419        """
420        Like fit, but more robust when each point has two or more attributes
421        """
422        try:
423            return self.fit(z)
424        except VectorShapeError, e:
425            # broadcasting is not supported.
426
427            #Convert input to Numeric arrays
428            z = array(z).astype(Float)
429           
430            #Build n x m interpolation matrix       
431            m = self.mesh.coordinates.shape[0] #Number of vertices
432            n = z.shape[1]               #Number of data points         
433
434            f = zeros((m,n), Float)
435            #f = sparse.dok_matrix() # even though it wont be sparse?
436           
437            for i in range(z.shape[1]):
438                f[:,i] = self.fit(z[:,i])
439            return f
440           
441       
442    def interpolate(self, f):
443        """Compute predicted values at data points implied in self.A.
444
445        The mesh values representing a smooth surface are
446        assumed to be specified in f. This argument could,
447        for example have been obtained from the method self.fit()
448       
449        Pre Condition:
450          self.A has been initialised
451       
452        Inputs:
453          f: Vector or array of data at the mesh vertices.
454          If f is an array, interpolation will be done for each column
455        """
456
457        try:
458            return self.A * f
459        except ValueError:
460            # We are here, so it probalby means that f is 2 dimensional
461            # and so will do each column separately due to problems in
462            # sparse matrix, 2d array multiplication
463            (N , M) = self.A.shape
464            f = array(f).astype(Float)
465            (m , n) = f.shape
466            #print "m",m
467            #print "M",M
468            if m != M :
469                #print "!!self.A",self.A
470                #print "!!self.A.todense() ",self.A.todense()
471                #print "self.get_A()",self.get_A()
472                #return dot(self.get_A(),f)
473                raise VectorShapeError, 'Mismatch between A and f dimensions'
474           
475            y = zeros( (N,n), Float)
476
477            for i in range(n):
478                y[:,i] = self.A * f[:,i]
479
480            #print "!!self.A.todense() ",self.A.todense()
481            return y
482     
483    #FIXME: We will need a method 'evaluate(self):' that will interpolate
484    #a computed surface living on the mesh onto a collection of
485    #arbitrary data points
486    #
487    #Precondition: self.fit(z) has stored its result in self.f.
488    #
489    #Input: data_points
490    #Algorithm:
491    #  1 Build a new temporary A matrix based on mesh and new data points
492    #  2 Apply it to self.f (return A*self.f)
493    #
494    # ON
495
496
497
498#-------------------------------------------------------------
499if __name__ == "__main__":
500    """
501    Load in a mesh and data points with attributes.
502    Fit the attributes to the mesh.
503    Save a new mesh file.
504    """
505    import os, sys
506    usage = "usage: %s mesh_input.tsh point.xya  mesh_output.tsh alpha" % os.path.basename(sys.argv[0])
507
508    if len(sys.argv) < 4:
509        print usage
510    else:
511        mesh_file = sys.argv[1]
512        point_file = sys.argv[2]
513        mesh_output_file = sys.argv[3]
514        if len(sys.argv) > 4:
515            alpha = sys.argv[4]
516        else:
517            alpha = DEFAULT_ALPHA
518        fit_to_mesh_file(mesh_file, point_file, mesh_output_file, alpha)
519       
Note: See TracBrowser for help on using the repository browser.