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

Last change on this file since 362 was 362, checked in by duncan, 20 years ago

adding attribute titles to loadASCII

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