[5897] | 1 | """Classes implementing general 2D triangular mesh with neighbour structure. |
---|
| 2 | |
---|
| 3 | This structure is purely geometrical. Anything relating to quantities |
---|
| 4 | or timestepping is implemented in subclass domain.py. |
---|
| 5 | |
---|
| 6 | Ole Nielsen, Stephen Roberts, Duncan Gray, Christopher Zoppou |
---|
| 7 | Geoscience Australia, 2004 |
---|
| 8 | """ |
---|
| 9 | |
---|
| 10 | from general_mesh import General_mesh |
---|
| 11 | from anuga.caching import cache |
---|
| 12 | from math import pi, sqrt |
---|
[6145] | 13 | |
---|
| 14 | import Numeric as num |
---|
[5897] | 15 | |
---|
| 16 | |
---|
| 17 | class Mesh(General_mesh): |
---|
| 18 | """Collection of triangular elements (purely geometric) |
---|
| 19 | |
---|
| 20 | A triangular element is defined in terms of three vertex ids, |
---|
| 21 | ordered counter clock-wise, |
---|
| 22 | each corresponding to a given coordinate set. |
---|
| 23 | Vertices from different elements can point to the same |
---|
| 24 | coordinate set. |
---|
| 25 | |
---|
| 26 | Coordinate sets are implemented as an N x 2 Numeric array containing |
---|
| 27 | x and y coordinates. |
---|
| 28 | |
---|
| 29 | |
---|
| 30 | To instantiate: |
---|
| 31 | Mesh(coordinates, triangles) |
---|
| 32 | |
---|
| 33 | where |
---|
| 34 | |
---|
| 35 | coordinates is either a list of 2-tuples or an Mx2 Numeric array of |
---|
| 36 | floats representing all x, y coordinates in the mesh. |
---|
| 37 | |
---|
| 38 | triangles is either a list of 3-tuples or an Nx3 Numeric array of |
---|
| 39 | integers representing indices of all vertices in the mesh. |
---|
| 40 | Each vertex is identified by its index i in [0, M-1]. |
---|
| 41 | |
---|
| 42 | |
---|
| 43 | Example: |
---|
| 44 | a = [0.0, 0.0] |
---|
| 45 | b = [0.0, 2.0] |
---|
| 46 | c = [2.0,0.0] |
---|
| 47 | e = [2.0, 2.0] |
---|
| 48 | |
---|
| 49 | points = [a, b, c, e] |
---|
| 50 | triangles = [ [1,0,2], [1,2,3] ] #bac, bce |
---|
| 51 | mesh = Mesh(points, triangles) |
---|
| 52 | |
---|
| 53 | #creates two triangles: bac and bce |
---|
| 54 | |
---|
| 55 | |
---|
| 56 | Mesh takes the optional third argument boundary which is a |
---|
| 57 | dictionary mapping from (element_id, edge_id) to boundary tag. |
---|
| 58 | The default value is None which will assign the default_boundary_tag |
---|
| 59 | as specified in config.py to all boundary edges. |
---|
| 60 | """ |
---|
| 61 | |
---|
| 62 | #FIXME: Maybe rename coordinates to points (as in a poly file) |
---|
| 63 | #But keep 'vertex_coordinates' |
---|
| 64 | |
---|
| 65 | #FIXME: Put in check for angles less than a set minimum |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | def __init__(self, coordinates, triangles, |
---|
| 69 | boundary=None, |
---|
| 70 | tagged_elements=None, |
---|
| 71 | geo_reference=None, |
---|
| 72 | number_of_full_nodes=None, |
---|
| 73 | number_of_full_triangles=None, |
---|
| 74 | use_inscribed_circle=False, |
---|
| 75 | verbose=False): |
---|
| 76 | """ |
---|
| 77 | Build triangles from x,y coordinates (sequence of 2-tuples or |
---|
| 78 | Mx2 Numeric array of floats) and triangles (sequence of 3-tuples |
---|
| 79 | or Nx3 Numeric array of non-negative integers). |
---|
| 80 | """ |
---|
| 81 | |
---|
| 82 | |
---|
| 83 | |
---|
| 84 | |
---|
| 85 | General_mesh.__init__(self, coordinates, triangles, |
---|
| 86 | number_of_full_nodes=\ |
---|
| 87 | number_of_full_nodes, |
---|
| 88 | number_of_full_triangles=\ |
---|
| 89 | number_of_full_triangles, |
---|
| 90 | geo_reference=geo_reference, |
---|
| 91 | verbose=verbose) |
---|
| 92 | |
---|
| 93 | if verbose: print 'Initialising mesh' |
---|
| 94 | |
---|
| 95 | N = len(self) #Number_of_triangles |
---|
| 96 | |
---|
| 97 | self.use_inscribed_circle = use_inscribed_circle |
---|
| 98 | |
---|
| 99 | #Allocate space for geometric quantities |
---|
[6145] | 100 | self.centroid_coordinates = num.zeros((N, 2), num.Float) |
---|
[5897] | 101 | |
---|
[6145] | 102 | self.radii = num.zeros(N, num.Float) |
---|
[5897] | 103 | |
---|
[6145] | 104 | self.neighbours = num.zeros((N, 3), num.Int) |
---|
| 105 | self.neighbour_edges = num.zeros((N, 3), num.Int) |
---|
| 106 | self.number_of_boundaries = num.zeros(N, num.Int) |
---|
| 107 | self.surrogate_neighbours = num.zeros((N, 3), num.Int) |
---|
[5897] | 108 | |
---|
| 109 | #Get x,y coordinates for all triangles and store |
---|
| 110 | V = self.vertex_coordinates # Relative coordinates |
---|
| 111 | |
---|
| 112 | #Initialise each triangle |
---|
| 113 | if verbose: print 'Mesh: Computing centroids and radii' |
---|
| 114 | for i in range(N): |
---|
| 115 | if verbose and i % ((N+10)/10) == 0: print '(%d/%d)' %(i, N) |
---|
| 116 | |
---|
| 117 | x0, y0 = V[3*i, :] |
---|
| 118 | x1, y1 = V[3*i+1, :] |
---|
| 119 | x2, y2 = V[3*i+2, :] |
---|
| 120 | |
---|
| 121 | #x0 = V[i, 0]; y0 = V[i, 1] |
---|
| 122 | #x1 = V[i, 2]; y1 = V[i, 3] |
---|
| 123 | #x2 = V[i, 4]; y2 = V[i, 5] |
---|
| 124 | |
---|
| 125 | #Compute centroid |
---|
[6174] | 126 | centroid = num.array([(x0 + x1 + x2)/3, (y0 + y1 + y2)/3], num.Float) |
---|
[5897] | 127 | self.centroid_coordinates[i] = centroid |
---|
| 128 | |
---|
| 129 | |
---|
| 130 | if self.use_inscribed_circle == False: |
---|
| 131 | #OLD code. Computed radii may exceed that of an |
---|
| 132 | #inscribed circle |
---|
| 133 | |
---|
| 134 | #Midpoints |
---|
[6174] | 135 | m0 = num.array([(x1 + x2)/2, (y1 + y2)/2], num.Float) |
---|
| 136 | m1 = num.array([(x0 + x2)/2, (y0 + y2)/2], num.Float) |
---|
| 137 | m2 = num.array([(x1 + x0)/2, (y1 + y0)/2], num.Float) |
---|
[5897] | 138 | |
---|
| 139 | #The radius is the distance from the centroid of |
---|
| 140 | #a triangle to the midpoint of the side of the triangle |
---|
| 141 | #closest to the centroid |
---|
[6145] | 142 | d0 = num.sqrt(num.sum( (centroid-m0)**2 )) |
---|
| 143 | d1 = num.sqrt(num.sum( (centroid-m1)**2 )) |
---|
| 144 | d2 = num.sqrt(num.sum( (centroid-m2)**2 )) |
---|
[5897] | 145 | |
---|
| 146 | self.radii[i] = min(d0, d1, d2) |
---|
| 147 | |
---|
| 148 | else: |
---|
| 149 | #NEW code added by Peter Row. True radius |
---|
| 150 | #of inscribed circle is computed |
---|
| 151 | |
---|
[6145] | 152 | a = num.sqrt((x0-x1)**2+(y0-y1)**2) |
---|
| 153 | b = num.sqrt((x1-x2)**2+(y1-y2)**2) |
---|
| 154 | c = num.sqrt((x2-x0)**2+(y2-y0)**2) |
---|
[5897] | 155 | |
---|
| 156 | self.radii[i]=2.0*self.areas[i]/(a+b+c) |
---|
| 157 | |
---|
| 158 | |
---|
| 159 | #Initialise Neighbours (-1 means that it is a boundary neighbour) |
---|
| 160 | self.neighbours[i, :] = [-1, -1, -1] |
---|
| 161 | |
---|
| 162 | #Initialise edge ids of neighbours |
---|
| 163 | #In case of boundaries this slot is not used |
---|
| 164 | self.neighbour_edges[i, :] = [-1, -1, -1] |
---|
| 165 | |
---|
| 166 | |
---|
| 167 | #Build neighbour structure |
---|
| 168 | if verbose: print 'Mesh: Building neigbour structure' |
---|
| 169 | self.build_neighbour_structure() |
---|
| 170 | |
---|
| 171 | #Build surrogate neighbour structure |
---|
| 172 | if verbose: print 'Mesh: Building surrogate neigbour structure' |
---|
| 173 | self.build_surrogate_neighbour_structure() |
---|
| 174 | |
---|
| 175 | #Build boundary dictionary mapping (id, edge) to symbolic tags |
---|
| 176 | if verbose: print 'Mesh: Building boundary dictionary' |
---|
| 177 | self.build_boundary_dictionary(boundary) |
---|
| 178 | |
---|
| 179 | #Build tagged element dictionary mapping (tag) to array of elements |
---|
| 180 | if verbose: print 'Mesh: Building tagged elements dictionary' |
---|
| 181 | self.build_tagged_elements_dictionary(tagged_elements) |
---|
| 182 | |
---|
| 183 | # Build a list of vertices that are not connected to any triangles |
---|
| 184 | self.lone_vertices = [] |
---|
| 185 | #Check that all vertices have been registered |
---|
| 186 | for node, count in enumerate(self.number_of_triangles_per_node): |
---|
| 187 | #msg = 'Node %d does not belong to an element.' %node |
---|
| 188 | #assert count > 0, msg |
---|
| 189 | if count == 0: |
---|
| 190 | self.lone_vertices.append(node) |
---|
| 191 | |
---|
| 192 | #Update boundary indices FIXME: OBSOLETE |
---|
| 193 | #self.build_boundary_structure() |
---|
| 194 | |
---|
| 195 | #FIXME check integrity? |
---|
| 196 | if verbose: print 'Mesh: Done' |
---|
| 197 | |
---|
| 198 | def __repr__(self): |
---|
| 199 | return General_mesh.__repr__(self) + ', %d boundary segments'\ |
---|
| 200 | %(len(self.boundary)) |
---|
| 201 | |
---|
| 202 | |
---|
| 203 | def set_to_inscribed_circle(self,safety_factor = 1): |
---|
| 204 | #FIXME phase out eventually |
---|
| 205 | N = self.number_of_triangles |
---|
| 206 | V = self.vertex_coordinates |
---|
| 207 | |
---|
| 208 | #initialising min and max ratio |
---|
| 209 | i=0 |
---|
| 210 | old_rad = self.radii[i] |
---|
| 211 | x0 = V[i, 0]; y0 = V[i, 1] |
---|
| 212 | x1 = V[i, 2]; y1 = V[i, 3] |
---|
| 213 | x2 = V[i, 4]; y2 = V[i, 5] |
---|
[6145] | 214 | a = num.sqrt((x0-x1)**2+(y0-y1)**2) |
---|
| 215 | b = num.sqrt((x1-x2)**2+(y1-y2)**2) |
---|
| 216 | c = num.sqrt((x2-x0)**2+(y2-y0)**2) |
---|
[5897] | 217 | ratio = old_rad/self.radii[i] |
---|
| 218 | max_ratio = ratio |
---|
| 219 | min_ratio = ratio |
---|
| 220 | |
---|
| 221 | for i in range(N): |
---|
| 222 | old_rad = self.radii[i] |
---|
| 223 | x0 = V[i, 0]; y0 = V[i, 1] |
---|
| 224 | x1 = V[i, 2]; y1 = V[i, 3] |
---|
| 225 | x2 = V[i, 4]; y2 = V[i, 5] |
---|
[6145] | 226 | a = num.sqrt((x0-x1)**2+(y0-y1)**2) |
---|
| 227 | b = num.sqrt((x1-x2)**2+(y1-y2)**2) |
---|
| 228 | c = num.sqrt((x2-x0)**2+(y2-y0)**2) |
---|
[5897] | 229 | self.radii[i]=self.areas[i]/(2*(a+b+c))*safety_factor |
---|
| 230 | ratio = old_rad/self.radii[i] |
---|
| 231 | if ratio >= max_ratio: max_ratio = ratio |
---|
| 232 | if ratio <= min_ratio: min_ratio = ratio |
---|
| 233 | return max_ratio,min_ratio |
---|
| 234 | |
---|
| 235 | |
---|
| 236 | |
---|
| 237 | def build_neighbour_structure(self): |
---|
| 238 | """Update all registered triangles to point to their neighbours. |
---|
| 239 | |
---|
| 240 | Also, keep a tally of the number of boundaries for each triangle |
---|
| 241 | |
---|
| 242 | Postconditions: |
---|
| 243 | neighbours and neighbour_edges is populated |
---|
| 244 | number_of_boundaries integer array is defined. |
---|
| 245 | """ |
---|
| 246 | |
---|
| 247 | #Step 1: |
---|
| 248 | #Build dictionary mapping from segments (2-tuple of points) |
---|
| 249 | #to left hand side edge (facing neighbouring triangle) |
---|
| 250 | |
---|
| 251 | N = len(self) #Number_of_triangles |
---|
| 252 | neighbourdict = {} |
---|
| 253 | for i in range(N): |
---|
| 254 | |
---|
| 255 | #Register all segments as keys mapping to current triangle |
---|
| 256 | #and segment id |
---|
| 257 | a = self.triangles[i, 0] |
---|
| 258 | b = self.triangles[i, 1] |
---|
| 259 | c = self.triangles[i, 2] |
---|
| 260 | if neighbourdict.has_key((a,b)): |
---|
| 261 | msg = "Edge 2 of triangle %d is duplicating edge %d of triangle %d.\n" %(i,neighbourdict[a,b][1],neighbourdict[a,b][0]) |
---|
| 262 | raise Exception, msg |
---|
| 263 | if neighbourdict.has_key((b,c)): |
---|
| 264 | msg = "Edge 0 of triangle %d is duplicating edge %d of triangle %d.\n" %(i,neighbourdict[b,c][1],neighbourdict[b,c][0]) |
---|
| 265 | raise msg |
---|
| 266 | if neighbourdict.has_key((c,a)): |
---|
| 267 | msg = "Edge 1 of triangle %d is duplicating edge %d of triangle %d.\n" %(i,neighbourdict[c,a][1],neighbourdict[c,a][0]) |
---|
| 268 | raise msg |
---|
| 269 | |
---|
| 270 | neighbourdict[a,b] = (i, 2) #(id, edge) |
---|
| 271 | neighbourdict[b,c] = (i, 0) #(id, edge) |
---|
| 272 | neighbourdict[c,a] = (i, 1) #(id, edge) |
---|
| 273 | |
---|
| 274 | |
---|
| 275 | #Step 2: |
---|
| 276 | #Go through triangles again, but this time |
---|
| 277 | #reverse direction of segments and lookup neighbours. |
---|
| 278 | for i in range(N): |
---|
| 279 | a = self.triangles[i, 0] |
---|
| 280 | b = self.triangles[i, 1] |
---|
| 281 | c = self.triangles[i, 2] |
---|
| 282 | |
---|
| 283 | self.number_of_boundaries[i] = 3 |
---|
| 284 | if neighbourdict.has_key((b,a)): |
---|
| 285 | self.neighbours[i, 2] = neighbourdict[b,a][0] |
---|
| 286 | self.neighbour_edges[i, 2] = neighbourdict[b,a][1] |
---|
| 287 | self.number_of_boundaries[i] -= 1 |
---|
| 288 | |
---|
| 289 | if neighbourdict.has_key((c,b)): |
---|
| 290 | self.neighbours[i, 0] = neighbourdict[c,b][0] |
---|
| 291 | self.neighbour_edges[i, 0] = neighbourdict[c,b][1] |
---|
| 292 | self.number_of_boundaries[i] -= 1 |
---|
| 293 | |
---|
| 294 | if neighbourdict.has_key((a,c)): |
---|
| 295 | self.neighbours[i, 1] = neighbourdict[a,c][0] |
---|
| 296 | self.neighbour_edges[i, 1] = neighbourdict[a,c][1] |
---|
| 297 | self.number_of_boundaries[i] -= 1 |
---|
| 298 | |
---|
| 299 | |
---|
| 300 | def build_surrogate_neighbour_structure(self): |
---|
| 301 | """Build structure where each triangle edge points to its neighbours |
---|
| 302 | if they exist. Otherwise point to the triangle itself. |
---|
| 303 | |
---|
| 304 | The surrogate neighbour structure is useful for computing gradients |
---|
| 305 | based on centroid values of neighbours. |
---|
| 306 | |
---|
| 307 | Precondition: Neighbour structure is defined |
---|
| 308 | Postcondition: |
---|
| 309 | Surrogate neighbour structure is defined: |
---|
| 310 | surrogate_neighbours: i0, i1, i2 where all i_k >= 0 point to |
---|
| 311 | triangles. |
---|
| 312 | |
---|
| 313 | """ |
---|
| 314 | |
---|
| 315 | N = len(self) #Number of triangles |
---|
| 316 | for i in range(N): |
---|
| 317 | #Find all neighbouring volumes that are not boundaries |
---|
| 318 | for k in range(3): |
---|
| 319 | if self.neighbours[i, k] < 0: |
---|
| 320 | self.surrogate_neighbours[i, k] = i #Point this triangle |
---|
| 321 | else: |
---|
| 322 | self.surrogate_neighbours[i, k] = self.neighbours[i, k] |
---|
| 323 | |
---|
| 324 | |
---|
| 325 | |
---|
| 326 | def build_boundary_dictionary(self, boundary = None): |
---|
| 327 | """Build or check the dictionary of boundary tags. |
---|
| 328 | self.boundary is a dictionary of tags, |
---|
| 329 | keyed by volume id and edge: |
---|
| 330 | { (id, edge): tag, ... } |
---|
| 331 | |
---|
| 332 | Postconditions: |
---|
| 333 | self.boundary is defined. |
---|
| 334 | """ |
---|
| 335 | |
---|
| 336 | from anuga.config import default_boundary_tag |
---|
| 337 | |
---|
| 338 | if boundary is None: |
---|
| 339 | boundary = {} |
---|
| 340 | for vol_id in range(len(self)): |
---|
| 341 | for edge_id in range(0, 3): |
---|
| 342 | if self.neighbours[vol_id, edge_id] < 0: |
---|
| 343 | boundary[(vol_id, edge_id)] = default_boundary_tag |
---|
| 344 | else: |
---|
| 345 | #Check that all keys in given boundary exist |
---|
| 346 | for vol_id, edge_id in boundary.keys(): |
---|
| 347 | msg = 'Segment (%d, %d) does not exist' %(vol_id, edge_id) |
---|
| 348 | a, b = self.neighbours.shape |
---|
| 349 | assert vol_id < a and edge_id < b, msg |
---|
| 350 | |
---|
| 351 | #FIXME: This assert violates internal boundaries (delete it) |
---|
| 352 | #msg = 'Segment (%d, %d) is not a boundary' %(vol_id, edge_id) |
---|
| 353 | #assert self.neighbours[vol_id, edge_id] < 0, msg |
---|
| 354 | |
---|
| 355 | #Check that all boundary segments are assigned a tag |
---|
| 356 | for vol_id in range(len(self)): |
---|
| 357 | for edge_id in range(0, 3): |
---|
| 358 | if self.neighbours[vol_id, edge_id] < 0: |
---|
| 359 | if not boundary.has_key( (vol_id, edge_id) ): |
---|
| 360 | msg = 'WARNING: Given boundary does not contain ' |
---|
| 361 | msg += 'tags for edge (%d, %d). '\ |
---|
| 362 | %(vol_id, edge_id) |
---|
| 363 | msg += 'Assigning default tag (%s).'\ |
---|
| 364 | %default_boundary_tag |
---|
| 365 | |
---|
| 366 | #FIXME: Print only as per verbosity |
---|
| 367 | #print msg |
---|
| 368 | |
---|
| 369 | #FIXME: Make this situation an error in the future |
---|
| 370 | #and make another function which will |
---|
| 371 | #enable default boundary-tags where |
---|
| 372 | #tags a not specified |
---|
| 373 | boundary[ (vol_id, edge_id) ] =\ |
---|
| 374 | default_boundary_tag |
---|
| 375 | |
---|
| 376 | |
---|
| 377 | |
---|
| 378 | self.boundary = boundary |
---|
| 379 | |
---|
| 380 | |
---|
| 381 | def build_tagged_elements_dictionary(self, tagged_elements = None): |
---|
| 382 | """Build the dictionary of element tags. |
---|
| 383 | self.tagged_elements is a dictionary of element arrays, |
---|
| 384 | keyed by tag: |
---|
| 385 | { (tag): [e1, e2, e3..] } |
---|
| 386 | |
---|
| 387 | Postconditions: |
---|
| 388 | self.element_tag is defined |
---|
| 389 | """ |
---|
| 390 | |
---|
| 391 | if tagged_elements is None: |
---|
| 392 | tagged_elements = {} |
---|
| 393 | else: |
---|
| 394 | #Check that all keys in given boundary exist |
---|
| 395 | for tag in tagged_elements.keys(): |
---|
[6166] | 396 | tagged_elements[tag] = num.array(tagged_elements[tag], num.Int) |
---|
[5897] | 397 | |
---|
| 398 | msg = 'Not all elements exist. ' |
---|
| 399 | assert max(tagged_elements[tag]) < len(self), msg |
---|
| 400 | #print "tagged_elements", tagged_elements |
---|
| 401 | self.tagged_elements = tagged_elements |
---|
| 402 | |
---|
| 403 | def build_boundary_structure(self): |
---|
| 404 | """Traverse boundary and |
---|
| 405 | enumerate neighbour indices from -1 and |
---|
| 406 | counting down. |
---|
| 407 | |
---|
| 408 | Precondition: |
---|
| 409 | self.boundary is defined. |
---|
| 410 | Post condition: |
---|
| 411 | neighbour array has unique negative indices for boundary |
---|
| 412 | boundary_segments array imposes an ordering on segments |
---|
| 413 | (not otherwise available from the dictionary) |
---|
| 414 | |
---|
| 415 | Note: If a segment is listed in the boundary dictionary |
---|
| 416 | it *will* become a boundary - even if there is a neighbouring triangle. |
---|
| 417 | This would be the case for internal boundaries |
---|
| 418 | """ |
---|
| 419 | |
---|
| 420 | #FIXME: Now Obsolete - maybe use some comments from here in |
---|
| 421 | #domain.set_boundary |
---|
| 422 | |
---|
| 423 | if self.boundary is None: |
---|
| 424 | msg = 'Boundary dictionary must be defined before ' |
---|
| 425 | msg += 'building boundary structure' |
---|
| 426 | raise msg |
---|
| 427 | |
---|
| 428 | |
---|
| 429 | self.boundary_segments = self.boundary.keys() |
---|
| 430 | self.boundary_segments.sort() |
---|
| 431 | |
---|
| 432 | index = -1 |
---|
| 433 | for id, edge in self.boundary_segments: |
---|
| 434 | |
---|
| 435 | #FIXME: One would detect internal boundaries as follows |
---|
| 436 | #if self.neighbours[id, edge] > -1: |
---|
| 437 | # print 'Internal boundary' |
---|
| 438 | |
---|
| 439 | self.neighbours[id, edge] = index |
---|
| 440 | index -= 1 |
---|
| 441 | |
---|
| 442 | |
---|
| 443 | def get_boundary_tags(self): |
---|
| 444 | """Return list of available boundary tags |
---|
| 445 | """ |
---|
| 446 | |
---|
| 447 | tags = {} |
---|
| 448 | for v in self.boundary.values(): |
---|
| 449 | tags[v] = 1 |
---|
| 450 | |
---|
| 451 | return tags.keys() |
---|
| 452 | |
---|
| 453 | |
---|
| 454 | def get_boundary_polygon(self, verbose=False): |
---|
| 455 | """Return bounding polygon for mesh (counter clockwise) |
---|
| 456 | |
---|
| 457 | Using the mesh boundary, derive a bounding polygon for this mesh. |
---|
| 458 | If multiple vertex values are present (vertices stored uniquely), |
---|
| 459 | the algorithm will select the path that contains the entire mesh. |
---|
| 460 | |
---|
| 461 | All points are in absolute UTM coordinates |
---|
| 462 | """ |
---|
| 463 | |
---|
| 464 | from anuga.utilities.numerical_tools import angle, ensure_numeric |
---|
| 465 | |
---|
| 466 | |
---|
| 467 | # Get mesh extent |
---|
| 468 | xmin, xmax, ymin, ymax = self.get_extent(absolute=True) |
---|
| 469 | pmin = ensure_numeric([xmin, ymin]) |
---|
| 470 | pmax = ensure_numeric([xmax, ymax]) |
---|
| 471 | |
---|
| 472 | |
---|
| 473 | # Assemble dictionary of boundary segments and choose starting point |
---|
| 474 | segments = {} |
---|
| 475 | inverse_segments = {} |
---|
| 476 | p0 = None |
---|
[6145] | 477 | mindist = num.sqrt(num.sum((pmax-pmin)**2)) # Start value across entire mesh |
---|
[5897] | 478 | for i, edge_id in self.boundary.keys(): |
---|
| 479 | # Find vertex ids for boundary segment |
---|
| 480 | if edge_id == 0: a = 1; b = 2 |
---|
| 481 | if edge_id == 1: a = 2; b = 0 |
---|
| 482 | if edge_id == 2: a = 0; b = 1 |
---|
| 483 | |
---|
| 484 | A = self.get_vertex_coordinate(i, a, absolute=True) # Start |
---|
| 485 | B = self.get_vertex_coordinate(i, b, absolute=True) # End |
---|
| 486 | |
---|
| 487 | |
---|
| 488 | # Take the point closest to pmin as starting point |
---|
| 489 | # Note: Could be arbitrary, but nice to have |
---|
| 490 | # a unique way of selecting |
---|
[6145] | 491 | dist_A = num.sqrt(num.sum((A-pmin)**2)) |
---|
| 492 | dist_B = num.sqrt(num.sum((B-pmin)**2)) |
---|
[5897] | 493 | |
---|
| 494 | # Find lower leftmost point |
---|
| 495 | if dist_A < mindist: |
---|
| 496 | mindist = dist_A |
---|
| 497 | p0 = A |
---|
| 498 | if dist_B < mindist: |
---|
| 499 | mindist = dist_B |
---|
| 500 | p0 = B |
---|
| 501 | |
---|
| 502 | |
---|
| 503 | # Sanity check |
---|
| 504 | if p0 is None: |
---|
| 505 | raise Exception('Impossible') |
---|
| 506 | |
---|
| 507 | |
---|
| 508 | # Register potential paths from A to B |
---|
| 509 | if not segments.has_key(tuple(A)): |
---|
| 510 | segments[tuple(A)] = [] # Empty list for candidate points |
---|
| 511 | |
---|
| 512 | segments[tuple(A)].append(B) |
---|
| 513 | |
---|
| 514 | |
---|
| 515 | # Start with smallest point and follow boundary (counter clock wise) |
---|
| 516 | polygon = [list(p0)]# Storage for final boundary polygon |
---|
| 517 | point_registry = {} # Keep track of storage to avoid multiple runs |
---|
| 518 | # around boundary. This will only be the case if |
---|
| 519 | # there are more than one candidate. |
---|
| 520 | # FIXME (Ole): Perhaps we can do away with polygon |
---|
| 521 | # and use only point_registry to save space. |
---|
| 522 | |
---|
| 523 | point_registry[tuple(p0)] = 0 |
---|
| 524 | |
---|
| 525 | while len(point_registry) < len(self.boundary): |
---|
| 526 | |
---|
| 527 | candidate_list = segments[tuple(p0)] |
---|
| 528 | if len(candidate_list) > 1: |
---|
| 529 | # Multiple points detected (this will be the case for meshes |
---|
| 530 | # with duplicate points as those used for discontinuous |
---|
| 531 | # triangles with vertices stored uniquely). |
---|
| 532 | # Take the candidate that is furthest to the clockwise |
---|
| 533 | # direction, as that will follow the boundary. |
---|
| 534 | # |
---|
| 535 | # This will also be the case for pathological triangles |
---|
| 536 | # that have no neighbours. |
---|
| 537 | |
---|
| 538 | if verbose: |
---|
| 539 | print 'Point %s has multiple candidates: %s'\ |
---|
| 540 | %(str(p0), candidate_list) |
---|
| 541 | |
---|
| 542 | # Check that previous are not in candidate list |
---|
| 543 | #for p in candidate_list: |
---|
| 544 | # assert not allclose(p0, p) |
---|
| 545 | |
---|
| 546 | # Choose vector against which all angles will be measured |
---|
| 547 | if len(polygon) > 1: |
---|
| 548 | v_prev = p0 - polygon[-2] # Vector that leads to p0 |
---|
| 549 | # from previous point |
---|
| 550 | else: |
---|
| 551 | # FIXME (Ole): What do we do if the first point has |
---|
| 552 | # multiple candidates? |
---|
| 553 | # Being the lower left corner, perhaps we can use the |
---|
| 554 | # vector [1, 0], but I really don't know if this is |
---|
| 555 | # completely watertight. |
---|
| 556 | v_prev = [1.0, 0.0] |
---|
| 557 | |
---|
| 558 | |
---|
| 559 | # Choose candidate with minimum angle |
---|
| 560 | minimum_angle = 2*pi |
---|
| 561 | for pc in candidate_list: |
---|
| 562 | |
---|
| 563 | vc = pc-p0 # Candidate vector (from p0 to candidate pt) |
---|
| 564 | |
---|
| 565 | # Angle between each candidate and the previous vector |
---|
| 566 | # in [-pi, pi] |
---|
| 567 | ac = angle(vc, v_prev) |
---|
| 568 | if ac > pi: |
---|
| 569 | # Give preference to angles on the right hand side |
---|
| 570 | # of v_prev |
---|
| 571 | # print 'pc = %s, changing angle from %f to %f'\ |
---|
| 572 | # %(pc, ac*180/pi, (ac-2*pi)*180/pi) |
---|
| 573 | ac = ac-2*pi |
---|
| 574 | |
---|
| 575 | # Take the minimal angle corresponding to the |
---|
| 576 | # rightmost vector |
---|
| 577 | if ac < minimum_angle: |
---|
| 578 | minimum_angle = ac |
---|
| 579 | p1 = pc # Best candidate |
---|
| 580 | |
---|
| 581 | |
---|
| 582 | if verbose is True: |
---|
| 583 | print ' Best candidate %s, angle %f'\ |
---|
| 584 | %(p1, minimum_angle*180/pi) |
---|
| 585 | |
---|
| 586 | else: |
---|
| 587 | p1 = candidate_list[0] |
---|
| 588 | |
---|
| 589 | |
---|
| 590 | if point_registry.has_key(tuple(p1)): |
---|
| 591 | # We have reached a point already visited. |
---|
| 592 | |
---|
[6145] | 593 | if num.allclose(p1, polygon[0]): |
---|
[5897] | 594 | # If it is the initial point, the polygon is complete. |
---|
| 595 | |
---|
| 596 | if verbose is True: |
---|
| 597 | print ' Stop criterion fulfilled at point %s' %p1 |
---|
| 598 | print polygon |
---|
| 599 | |
---|
| 600 | # We have completed the boundary polygon - yeehaa |
---|
| 601 | break |
---|
| 602 | else: |
---|
| 603 | # The point already visited is not the initial point |
---|
| 604 | # This would be a pathological triangle, but the |
---|
| 605 | # algorithm must be able to deal with this |
---|
| 606 | pass |
---|
| 607 | |
---|
| 608 | else: |
---|
| 609 | # We are still finding new points on the boundary |
---|
| 610 | point_registry[tuple(p1)] = len(point_registry) |
---|
| 611 | |
---|
| 612 | polygon.append(list(p1)) # De-Numeric each point :-) |
---|
| 613 | p0 = p1 |
---|
| 614 | |
---|
| 615 | |
---|
| 616 | return polygon |
---|
| 617 | |
---|
| 618 | |
---|
| 619 | def check_integrity(self): |
---|
| 620 | """Check that triangles are internally consistent e.g. |
---|
| 621 | that area corresponds to edgelengths, that vertices |
---|
| 622 | are arranged in a counter-clockwise order, etc etc |
---|
| 623 | Neighbour structure will be checked by class Mesh |
---|
| 624 | """ |
---|
| 625 | |
---|
| 626 | from anuga.config import epsilon |
---|
| 627 | from anuga.utilities.numerical_tools import anglediff |
---|
| 628 | |
---|
| 629 | N = len(self) |
---|
| 630 | |
---|
| 631 | # Get x,y coordinates for all vertices for all triangles |
---|
| 632 | V = self.get_vertex_coordinates() |
---|
| 633 | |
---|
| 634 | # Check each triangle |
---|
| 635 | for i in range(N): |
---|
| 636 | |
---|
| 637 | x0, y0 = V[3*i, :] |
---|
| 638 | x1, y1 = V[3*i+1, :] |
---|
| 639 | x2, y2 = V[3*i+2, :] |
---|
| 640 | |
---|
| 641 | # Check that area hasn't been compromised |
---|
| 642 | area = self.areas[i] |
---|
| 643 | ref = abs((x1*y0-x0*y1)+(x2*y1-x1*y2)+(x0*y2-x2*y0))/2 |
---|
| 644 | msg = 'Wrong area for vertex coordinates: %f %f %f %f %f %f'\ |
---|
| 645 | %(x0,y0,x1,y1,x2,y2) |
---|
| 646 | assert abs((area - ref)/area) < epsilon, msg |
---|
| 647 | |
---|
| 648 | # Check that points are arranged in counter clock-wise order |
---|
| 649 | v0 = [x1-x0, y1-y0] |
---|
| 650 | v1 = [x2-x1, y2-y1] |
---|
| 651 | v2 = [x0-x2, y0-y2] |
---|
| 652 | a0 = anglediff(v1, v0) |
---|
| 653 | a1 = anglediff(v2, v1) |
---|
| 654 | a2 = anglediff(v0, v2) |
---|
| 655 | |
---|
| 656 | msg = '''Vertices (%s,%s), (%s,%s), (%s,%s) are not arranged |
---|
| 657 | in counter clockwise order''' %(x0, y0, x1, y1, x2, y2) |
---|
| 658 | assert a0 < pi and a1 < pi and a2 < pi, msg |
---|
| 659 | |
---|
| 660 | # Check that normals are orthogonal to edge vectors |
---|
| 661 | # Note that normal[k] lies opposite vertex k |
---|
| 662 | |
---|
| 663 | normal0 = self.normals[i, 0:2] |
---|
| 664 | normal1 = self.normals[i, 2:4] |
---|
| 665 | normal2 = self.normals[i, 4:6] |
---|
| 666 | |
---|
| 667 | for u, v in [ (v0, normal2), (v1, normal0), (v2, normal1) ]: |
---|
| 668 | |
---|
| 669 | # Normalise |
---|
[6145] | 670 | l_u = num.sqrt(u[0]*u[0] + u[1]*u[1]) |
---|
| 671 | l_v = num.sqrt(v[0]*v[0] + v[1]*v[1]) |
---|
[5897] | 672 | |
---|
| 673 | msg = 'Normal vector in triangle %d does not have unit length' %i |
---|
[6145] | 674 | assert num.allclose(l_v, 1), msg |
---|
[5897] | 675 | |
---|
| 676 | x = (u[0]*v[0] + u[1]*v[1])/l_u # Inner product |
---|
| 677 | |
---|
| 678 | msg = 'Normal vector (%f,%f) is not perpendicular to' %tuple(v) |
---|
| 679 | msg += ' edge (%f,%f) in triangle %d.' %(tuple(u) + (i,)) |
---|
| 680 | msg += ' Inner product is %e.' %x |
---|
| 681 | assert x < epsilon, msg |
---|
| 682 | |
---|
| 683 | |
---|
| 684 | |
---|
| 685 | |
---|
| 686 | #Check neighbour structure |
---|
| 687 | for i in range(N): |
---|
| 688 | # For each triangle |
---|
| 689 | |
---|
| 690 | for k, neighbour_id in enumerate(self.neighbours[i,:]): |
---|
| 691 | |
---|
| 692 | #Assert that my neighbour's neighbour is me |
---|
| 693 | #Boundaries need not fulfill this |
---|
| 694 | if neighbour_id >= 0: |
---|
| 695 | edge = self.neighbour_edges[i, k] |
---|
| 696 | msg = 'Triangle %d has neighbour %d but it does not point back. \n' %(i,neighbour_id) |
---|
| 697 | msg += 'Only points to (%s)' %(self.neighbours[neighbour_id,:]) |
---|
| 698 | assert self.neighbours[neighbour_id, edge] == i ,msg |
---|
| 699 | |
---|
| 700 | |
---|
| 701 | |
---|
| 702 | #Check that all boundaries have |
---|
| 703 | # unique, consecutive, negative indices |
---|
| 704 | |
---|
| 705 | #L = len(self.boundary) |
---|
| 706 | #for i in range(L): |
---|
| 707 | # id, edge = self.boundary_segments[i] |
---|
| 708 | # assert self.neighbours[id, edge] == -i-1 |
---|
| 709 | |
---|
| 710 | |
---|
| 711 | #NOTE: This assert doesn't hold true if there are internal boundaries |
---|
| 712 | #FIXME: Look into this further. |
---|
| 713 | #FIXME (Ole): In pyvolution mark 3 this is OK again |
---|
| 714 | #NOTE: No longer works because neighbour structure is modified by |
---|
| 715 | # domain set_boundary. |
---|
| 716 | #for id, edge in self.boundary: |
---|
| 717 | # assert self.neighbours[id,edge] < 0 |
---|
| 718 | # |
---|
| 719 | #NOTE (Ole): I reckon this was resolved late 2004? |
---|
| 720 | # |
---|
| 721 | #See domain.set_boundary |
---|
| 722 | |
---|
| 723 | |
---|
| 724 | |
---|
| 725 | #Check integrity of inverted triangle structure |
---|
| 726 | |
---|
| 727 | V = self.vertex_value_indices[:] #Take a copy |
---|
[6145] | 728 | V = num.sort(V) |
---|
| 729 | assert num.allclose(V, range(3*N)) |
---|
[5897] | 730 | |
---|
[6145] | 731 | assert num.sum(self.number_of_triangles_per_node) ==\ |
---|
| 732 | len(self.vertex_value_indices) |
---|
[5897] | 733 | |
---|
| 734 | # Check number of triangles per node |
---|
| 735 | count = [0]*self.number_of_nodes |
---|
| 736 | for triangle in self.triangles: |
---|
| 737 | for i in triangle: |
---|
| 738 | count[i] += 1 |
---|
| 739 | |
---|
[6145] | 740 | assert num.allclose(count, self.number_of_triangles_per_node) |
---|
[5897] | 741 | |
---|
| 742 | |
---|
| 743 | # Check integrity of vertex_value_indices |
---|
| 744 | current_node = 0 |
---|
| 745 | k = 0 # Track triangles touching on node |
---|
| 746 | for index in self.vertex_value_indices: |
---|
| 747 | |
---|
| 748 | if self.number_of_triangles_per_node[current_node] == 0: |
---|
| 749 | # Node is lone - i.e. not part of the mesh |
---|
| 750 | continue |
---|
| 751 | |
---|
| 752 | k += 1 |
---|
| 753 | |
---|
| 754 | volume_id = index / 3 |
---|
| 755 | vertex_id = index % 3 |
---|
| 756 | |
---|
| 757 | msg = 'Triangle %d, vertex %d points to %d. Should have been %d'\ |
---|
| 758 | %(volume_id, vertex_id, self.triangles[volume_id, vertex_id], current_node) |
---|
| 759 | assert self.triangles[volume_id, vertex_id] == current_node, msg |
---|
| 760 | |
---|
| 761 | if self.number_of_triangles_per_node[current_node] == k: |
---|
| 762 | # Move on to next node |
---|
| 763 | k = 0 |
---|
| 764 | current_node += 1 |
---|
| 765 | |
---|
| 766 | |
---|
| 767 | def get_lone_vertices(self): |
---|
| 768 | """Return a list of vertices that are not connected to any triangles. |
---|
| 769 | |
---|
| 770 | """ |
---|
| 771 | return self.lone_vertices |
---|
| 772 | |
---|
| 773 | def get_centroid_coordinates(self, absolute=False): |
---|
| 774 | """Return all centroid coordinates. |
---|
| 775 | Return all centroid coordinates for all triangles as an Nx2 array |
---|
| 776 | (ordered as x0, y0 for each triangle) |
---|
| 777 | |
---|
| 778 | Boolean keyword argument absolute determines whether coordinates |
---|
| 779 | are to be made absolute by taking georeference into account |
---|
| 780 | Default is False as many parts of ANUGA expects relative coordinates. |
---|
| 781 | """ |
---|
| 782 | |
---|
| 783 | V = self.centroid_coordinates |
---|
| 784 | if absolute is True: |
---|
| 785 | if not self.geo_reference.is_absolute(): |
---|
| 786 | V = self.geo_reference.get_absolute(V) |
---|
| 787 | |
---|
| 788 | return V |
---|
| 789 | |
---|
| 790 | |
---|
| 791 | def get_radii(self): |
---|
| 792 | """Return all radii. |
---|
| 793 | Return radius of inscribed cirle for all triangles |
---|
| 794 | """ |
---|
| 795 | return self.radii |
---|
| 796 | |
---|
| 797 | |
---|
| 798 | |
---|
| 799 | def statistics(self): |
---|
| 800 | """Output statistics about mesh |
---|
| 801 | """ |
---|
| 802 | |
---|
| 803 | from anuga.utilities.numerical_tools import histogram, create_bins |
---|
| 804 | |
---|
| 805 | vertex_coordinates = self.vertex_coordinates # Relative coordinates |
---|
| 806 | areas = self.areas |
---|
| 807 | x = vertex_coordinates[:,0] |
---|
| 808 | y = vertex_coordinates[:,1] |
---|
| 809 | |
---|
| 810 | |
---|
| 811 | #Setup 10 bins for area histogram |
---|
| 812 | bins = create_bins(areas, 10) |
---|
| 813 | #m = max(areas) |
---|
| 814 | #bins = arange(0., m, m/10) |
---|
| 815 | hist = histogram(areas, bins) |
---|
| 816 | |
---|
| 817 | str = '------------------------------------------------\n' |
---|
| 818 | str += 'Mesh statistics:\n' |
---|
| 819 | str += ' Number of triangles = %d\n' %len(self) |
---|
| 820 | str += ' Extent [m]:\n' |
---|
| 821 | str += ' x in [%f, %f]\n' %(min(x), max(x)) |
---|
| 822 | str += ' y in [%f, %f]\n' %(min(y), max(y)) |
---|
| 823 | str += ' Areas [m^2]:\n' |
---|
| 824 | str += ' A in [%f, %f]\n' %(min(areas), max(areas)) |
---|
| 825 | str += ' number of distinct areas: %d\n' %(len(areas)) |
---|
| 826 | str += ' Histogram:\n' |
---|
| 827 | |
---|
| 828 | hi = bins[0] |
---|
| 829 | for i, count in enumerate(hist): |
---|
| 830 | lo = hi |
---|
| 831 | if i+1 < len(bins): |
---|
| 832 | #Open upper interval |
---|
| 833 | hi = bins[i+1] |
---|
| 834 | str += ' [%f, %f[: %d\n' %(lo, hi, count) |
---|
| 835 | else: |
---|
| 836 | #Closed upper interval |
---|
| 837 | hi = max(areas) |
---|
| 838 | str += ' [%f, %f]: %d\n' %(lo, hi, count) |
---|
| 839 | |
---|
| 840 | N = len(areas) |
---|
| 841 | if N > 10: |
---|
| 842 | str += ' Percentiles (10%):\n' |
---|
| 843 | areas = areas.tolist() |
---|
| 844 | areas.sort() |
---|
| 845 | |
---|
| 846 | k = 0 |
---|
| 847 | lower = min(areas) |
---|
| 848 | for i, a in enumerate(areas): |
---|
| 849 | if i % (N/10) == 0 and i != 0: #For every 10% of the sorted areas |
---|
| 850 | str += ' %d triangles in [%f, %f]\n' %(i-k, lower, a) |
---|
| 851 | lower = a |
---|
| 852 | k = i |
---|
| 853 | |
---|
| 854 | str += ' %d triangles in [%f, %f]\n'\ |
---|
| 855 | %(N-k, lower, max(areas)) |
---|
| 856 | |
---|
| 857 | |
---|
| 858 | str += ' Boundary:\n' |
---|
| 859 | str += ' Number of boundary segments == %d\n' %(len(self.boundary)) |
---|
| 860 | str += ' Boundary tags == %s\n' %self.get_boundary_tags() |
---|
| 861 | str += '------------------------------------------------\n' |
---|
| 862 | |
---|
| 863 | |
---|
| 864 | return str |
---|
| 865 | |
---|
| 866 | |
---|
| 867 | def get_triangle_containing_point(self, point): |
---|
| 868 | """Return triangle id for triangle containing specifiend point (x,y) |
---|
| 869 | |
---|
| 870 | If point isn't within mesh, raise exception |
---|
| 871 | |
---|
| 872 | """ |
---|
| 873 | |
---|
| 874 | # FIXME(Ole): This function is currently brute force |
---|
| 875 | # because I needed it for diagnostics. |
---|
| 876 | # We should make it fast - probably based on the |
---|
| 877 | # quad tree structure. |
---|
| 878 | from anuga.utilities.polygon import is_outside_polygon,\ |
---|
| 879 | is_inside_polygon |
---|
| 880 | |
---|
| 881 | polygon = self.get_boundary_polygon() |
---|
| 882 | #print polygon, point |
---|
| 883 | |
---|
| 884 | if is_outside_polygon(point, polygon): |
---|
| 885 | msg = 'Point %s is outside mesh' %str(point) |
---|
| 886 | raise Exception, msg |
---|
| 887 | |
---|
| 888 | |
---|
| 889 | V = self.get_vertex_coordinates(absolute=True) |
---|
| 890 | |
---|
| 891 | # FIXME: Horrible brute force |
---|
| 892 | for i, triangle in enumerate(self.triangles): |
---|
| 893 | poly = V[3*i:3*i+3] |
---|
| 894 | #print i, poly |
---|
| 895 | |
---|
| 896 | if is_inside_polygon(point, poly, closed=True): |
---|
| 897 | return i |
---|
| 898 | |
---|
| 899 | return |
---|
| 900 | |
---|
| 901 | |
---|
| 902 | |
---|
| 903 | def get_intersecting_segments(self, polyline, |
---|
| 904 | use_cache=False, |
---|
| 905 | verbose=False): |
---|
| 906 | """Find edges intersected by polyline |
---|
| 907 | |
---|
| 908 | Input: |
---|
| 909 | polyline - list of points forming a segmented line |
---|
| 910 | use_cache |
---|
| 911 | verbose |
---|
| 912 | |
---|
| 913 | Output: |
---|
| 914 | list of instances of class Triangle_intersection |
---|
| 915 | |
---|
| 916 | The polyline may break inside any triangle causing multiple |
---|
| 917 | segments per triangle - consequently the same triangle may |
---|
| 918 | appear in several entries. |
---|
| 919 | |
---|
| 920 | If a polyline segment coincides with a triangle edge, |
---|
| 921 | the the entire shared segment will be used. |
---|
| 922 | Onle one of the triangles thus intersected will be used and that |
---|
| 923 | is the first one encountered. |
---|
| 924 | |
---|
| 925 | Intersections with single vertices are ignored. |
---|
| 926 | |
---|
| 927 | Resulting segments are unsorted |
---|
| 928 | """ |
---|
| 929 | |
---|
| 930 | V = self.get_vertex_coordinates() |
---|
| 931 | N = len(self) |
---|
| 932 | |
---|
| 933 | # Adjust polyline to mesh spatial origin |
---|
| 934 | polyline = self.geo_reference.get_relative(polyline) |
---|
| 935 | |
---|
| 936 | if use_cache is True: |
---|
| 937 | segments = cache(get_intersecting_segments, |
---|
| 938 | (V, N, polyline), |
---|
| 939 | {'verbose': verbose}, |
---|
| 940 | verbose=verbose) |
---|
| 941 | else: |
---|
| 942 | segments = get_intersecting_segments(V, N, polyline, |
---|
| 943 | verbose=verbose) |
---|
| 944 | |
---|
| 945 | |
---|
| 946 | return segments |
---|
| 947 | |
---|
| 948 | |
---|
| 949 | |
---|
| 950 | def get_triangle_neighbours(self, tri_id): |
---|
| 951 | """ Given a triangle id, Return an array of the |
---|
| 952 | 3 neighbour triangle id's. |
---|
| 953 | |
---|
| 954 | Negative returned triangle id's represent a boundary as a neighbour. |
---|
| 955 | |
---|
| 956 | If the given triangle id is bad, return an empty list. |
---|
| 957 | """ |
---|
| 958 | |
---|
| 959 | try: |
---|
| 960 | return self.neighbours[tri_id,:] |
---|
| 961 | except IndexError: |
---|
| 962 | return [] |
---|
| 963 | |
---|
| 964 | |
---|
| 965 | def get_interpolation_object(self): |
---|
| 966 | """Get object I that will allow linear interpolation using this mesh |
---|
| 967 | |
---|
| 968 | This is a time consuming process but it needs only to be |
---|
| 969 | once for the mesh. |
---|
| 970 | |
---|
| 971 | Interpolation can then be done using |
---|
| 972 | |
---|
| 973 | result = I.interpolate_block(vertex_values, interpolation_points) |
---|
| 974 | |
---|
| 975 | where vertex values have been obtained from a quantity using |
---|
| 976 | vertex_values, triangles = self.get_vertex_values() |
---|
| 977 | """ |
---|
| 978 | |
---|
| 979 | if hasattr(self, 'interpolation_object'): |
---|
| 980 | I = self.interpolation_object |
---|
| 981 | else: |
---|
| 982 | from anuga.fit_interpolate.interpolate import Interpolate |
---|
| 983 | |
---|
| 984 | # Get discontinuous mesh - this will match internal |
---|
| 985 | # representation of vertex values |
---|
| 986 | triangles = self.get_disconnected_triangles() |
---|
| 987 | vertex_coordinates = self.get_vertex_coordinates() |
---|
| 988 | |
---|
| 989 | I = Interpolate(vertex_coordinates, triangles) |
---|
| 990 | self.interpolation_object = I |
---|
| 991 | |
---|
| 992 | return I |
---|
| 993 | |
---|
| 994 | |
---|
| 995 | class Triangle_intersection: |
---|
| 996 | """Store information about line segments intersecting a triangle |
---|
| 997 | |
---|
| 998 | Attributes are |
---|
| 999 | |
---|
| 1000 | segment: Line segment intersecting triangle [[x0,y0], [x1, y1]] |
---|
| 1001 | normal: [a,b] right hand normal to segment |
---|
| 1002 | length: Length of intersecting segment |
---|
| 1003 | triangle_id: id (in mesh) of triangle being intersected |
---|
| 1004 | |
---|
| 1005 | """ |
---|
| 1006 | |
---|
| 1007 | |
---|
| 1008 | def __init__(self, |
---|
| 1009 | segment=None, |
---|
| 1010 | normal=None, |
---|
| 1011 | length=None, |
---|
| 1012 | triangle_id=None): |
---|
| 1013 | self.segment = segment |
---|
| 1014 | self.normal = normal |
---|
| 1015 | self.length = length |
---|
| 1016 | self.triangle_id = triangle_id |
---|
| 1017 | |
---|
| 1018 | |
---|
| 1019 | def __repr__(self): |
---|
| 1020 | s = 'Triangle_intersection(' |
---|
| 1021 | s += 'segment=%s, normal=%s, length=%s, triangle_id=%s)'\ |
---|
| 1022 | %(self.segment, |
---|
| 1023 | self.normal, |
---|
| 1024 | self.length, |
---|
| 1025 | self.triangle_id) |
---|
| 1026 | |
---|
| 1027 | return s |
---|
| 1028 | |
---|
| 1029 | |
---|
| 1030 | |
---|
| 1031 | def _get_intersecting_segments(V, N, line, |
---|
| 1032 | verbose=False): |
---|
| 1033 | """Find edges intersected by line |
---|
| 1034 | |
---|
| 1035 | Input: |
---|
| 1036 | V: Vertex coordinates as obtained by mesh.get_vertex_coordinates() |
---|
| 1037 | N: Number of triangles in mesh |
---|
| 1038 | line - list of two points forming a segmented line |
---|
| 1039 | verbose |
---|
| 1040 | Output: |
---|
| 1041 | list of instances of class Triangle_intersection |
---|
| 1042 | |
---|
| 1043 | This method is used by the public method |
---|
| 1044 | get_intersecting_segments(self, polyline) which also contains |
---|
| 1045 | more documentation. |
---|
| 1046 | """ |
---|
| 1047 | |
---|
| 1048 | from anuga.utilities.polygon import intersection |
---|
| 1049 | from anuga.utilities.polygon import is_inside_polygon |
---|
| 1050 | |
---|
| 1051 | msg = 'Line segment must contain exactly two points' |
---|
| 1052 | assert len(line) == 2, msg |
---|
| 1053 | |
---|
| 1054 | # Origin of intersecting line to be used for |
---|
| 1055 | # establishing direction |
---|
| 1056 | xi0 = line[0][0] |
---|
| 1057 | eta0 = line[0][1] |
---|
| 1058 | |
---|
| 1059 | |
---|
| 1060 | # Check intersection with edge segments for all triangles |
---|
| 1061 | # FIXME (Ole): This should be implemented in C |
---|
| 1062 | triangle_intersections={} # Keep track of segments already done |
---|
| 1063 | for i in range(N): |
---|
| 1064 | # Get nodes and edge segments for each triangle |
---|
| 1065 | x0, y0 = V[3*i, :] |
---|
| 1066 | x1, y1 = V[3*i+1, :] |
---|
| 1067 | x2, y2 = V[3*i+2, :] |
---|
| 1068 | |
---|
| 1069 | |
---|
| 1070 | edge_segments = [[[x0,y0], [x1, y1]], |
---|
| 1071 | [[x1,y1], [x2, y2]], |
---|
| 1072 | [[x2,y2], [x0, y0]]] |
---|
| 1073 | |
---|
| 1074 | # Find segments that are intersected by line |
---|
| 1075 | |
---|
| 1076 | intersections = {} # Use dictionary to record points only once |
---|
| 1077 | for edge in edge_segments: |
---|
| 1078 | |
---|
| 1079 | status, value = intersection(line, edge) |
---|
| 1080 | #if value is not None: print 'Triangle %d, Got' %i, status, value |
---|
| 1081 | |
---|
| 1082 | if status == 1: |
---|
| 1083 | # Normal intersection of one edge or vertex |
---|
| 1084 | intersections[tuple(value)] = i |
---|
| 1085 | |
---|
| 1086 | # Exclude singular intersections with vertices |
---|
| 1087 | #if not(allclose(value, edge[0]) or\ |
---|
| 1088 | # allclose(value, edge[1])): |
---|
| 1089 | # intersections.append(value) |
---|
| 1090 | |
---|
| 1091 | if status == 2: |
---|
| 1092 | # Edge is sharing a segment with line |
---|
| 1093 | |
---|
| 1094 | # This is usually covered by the two |
---|
| 1095 | # vertices that would have been picked up |
---|
| 1096 | # under status == 1. |
---|
| 1097 | # However, if coinciding line stops partway |
---|
| 1098 | # along this edge, it will be recorded here. |
---|
| 1099 | intersections[tuple(value[0,:])] = i |
---|
| 1100 | intersections[tuple(value[1,:])] = i |
---|
| 1101 | |
---|
| 1102 | |
---|
| 1103 | if len(intersections) == 1: |
---|
| 1104 | # Check if either line end point lies fully within this triangle |
---|
| 1105 | # If this is the case accept that as one end of the intersecting |
---|
| 1106 | # segment |
---|
| 1107 | |
---|
| 1108 | poly = V[3*i:3*i+3] |
---|
| 1109 | if is_inside_polygon(line[1], poly, closed=False): |
---|
| 1110 | intersections[tuple(line[1])] = i |
---|
| 1111 | elif is_inside_polygon(line[0], poly, closed=False): |
---|
| 1112 | intersections[tuple(line[0])] = i |
---|
| 1113 | else: |
---|
| 1114 | # Ignore situations where one vertex is touch, for instance |
---|
| 1115 | continue |
---|
| 1116 | |
---|
| 1117 | |
---|
| 1118 | msg = 'There can be only two or no intersections' |
---|
| 1119 | assert len(intersections) in [0,2], msg |
---|
| 1120 | |
---|
| 1121 | |
---|
| 1122 | if len(intersections) == 2: |
---|
| 1123 | |
---|
| 1124 | # Calculate attributes for this segment |
---|
| 1125 | |
---|
| 1126 | |
---|
| 1127 | # End points of intersecting segment |
---|
| 1128 | points = intersections.keys() |
---|
| 1129 | x0, y0 = points[0] |
---|
| 1130 | x1, y1 = points[1] |
---|
| 1131 | |
---|
| 1132 | |
---|
| 1133 | # Determine which end point is closer to the origin of the line |
---|
| 1134 | # This is necessary for determining the direction of |
---|
| 1135 | # the line and the normals |
---|
| 1136 | |
---|
| 1137 | # Distances from line origin to the two intersections |
---|
[6174] | 1138 | z0 = num.array([x0 - xi0, y0 - eta0], num.Float) |
---|
| 1139 | z1 = num.array([x1 - xi0, y1 - eta0], num.Float) |
---|
[6145] | 1140 | d0 = num.sqrt(num.sum(z0**2)) |
---|
| 1141 | d1 = num.sqrt(num.sum(z1**2)) |
---|
[5897] | 1142 | |
---|
| 1143 | if d1 < d0: |
---|
| 1144 | # Swap |
---|
| 1145 | xi, eta = x0, y0 |
---|
| 1146 | x0, y0 = x1, y1 |
---|
| 1147 | x1, y1 = xi, eta |
---|
| 1148 | |
---|
| 1149 | # (x0,y0) is now the origin of the intersecting segment |
---|
| 1150 | |
---|
| 1151 | |
---|
| 1152 | # Normal direction: |
---|
| 1153 | # Right hand side relative to line direction |
---|
[6174] | 1154 | vector = num.array([x1 - x0, y1 - y0], num.Float) # Segment vector |
---|
[6145] | 1155 | length = num.sqrt(num.sum(vector**2)) # Segment length |
---|
[6174] | 1156 | normal = num.array([vector[1], -vector[0]], num.Float)/length |
---|
[5897] | 1157 | |
---|
| 1158 | |
---|
| 1159 | segment = ((x0,y0), (x1, y1)) |
---|
| 1160 | T = Triangle_intersection(segment=segment, |
---|
| 1161 | normal=normal, |
---|
| 1162 | length=length, |
---|
| 1163 | triangle_id=i) |
---|
| 1164 | |
---|
| 1165 | |
---|
| 1166 | # Add segment unless it was done earlier |
---|
| 1167 | if not triangle_intersections.has_key(segment): |
---|
| 1168 | triangle_intersections[segment] = T |
---|
| 1169 | |
---|
| 1170 | |
---|
| 1171 | # Return segments as a list |
---|
| 1172 | return triangle_intersections.values() |
---|
| 1173 | |
---|
| 1174 | |
---|
| 1175 | def get_intersecting_segments(V, N, polyline, |
---|
| 1176 | verbose=False): |
---|
| 1177 | """Internal function to find edges intersected by Polyline |
---|
| 1178 | |
---|
| 1179 | Input: |
---|
| 1180 | V: Vertex coordinates as obtained by mesh.get_vertex_coordinates() |
---|
| 1181 | N: Number of triangles in mesh |
---|
| 1182 | polyline - list of points forming a segmented line |
---|
| 1183 | verbose |
---|
| 1184 | Output: |
---|
| 1185 | list of instances of class Triangle_intersection |
---|
| 1186 | |
---|
| 1187 | This method is used by the public method |
---|
| 1188 | get_intersecting_segments(self, polyline) which also contains |
---|
| 1189 | more documentation. |
---|
| 1190 | """ |
---|
| 1191 | |
---|
| 1192 | msg = 'Polyline must contain at least two points' |
---|
| 1193 | assert len(polyline) >= 2, msg |
---|
| 1194 | |
---|
| 1195 | |
---|
| 1196 | # For all segments in polyline |
---|
| 1197 | triangle_intersections = [] |
---|
| 1198 | for i, point0 in enumerate(polyline[:-1]): |
---|
| 1199 | |
---|
| 1200 | point1 = polyline[i+1] |
---|
| 1201 | if verbose: |
---|
| 1202 | print 'Extracting mesh intersections from line:', |
---|
| 1203 | print '(%.2f, %.2f) - (%.2f, %.2f)' %(point0[0], point0[1], |
---|
| 1204 | point1[0], point1[1]) |
---|
| 1205 | |
---|
| 1206 | line = [point0, point1] |
---|
| 1207 | triangle_intersections += _get_intersecting_segments(V, N, line, |
---|
| 1208 | verbose=verbose) |
---|
| 1209 | |
---|
| 1210 | |
---|
| 1211 | msg = 'No segments found' |
---|
| 1212 | assert len(triangle_intersections) > 0, msg |
---|
| 1213 | |
---|
| 1214 | |
---|
| 1215 | return triangle_intersections |
---|
| 1216 | |
---|
| 1217 | |
---|
| 1218 | |
---|
| 1219 | |
---|
| 1220 | |
---|
| 1221 | def segment_midpoints(segments): |
---|
| 1222 | """Calculate midpoints of all segments |
---|
| 1223 | |
---|
| 1224 | Inputs: |
---|
| 1225 | segments: List of instances of class Segment |
---|
| 1226 | |
---|
| 1227 | Ouputs: |
---|
| 1228 | midpoints: List of points |
---|
| 1229 | """ |
---|
| 1230 | |
---|
| 1231 | midpoints = [] |
---|
| 1232 | msg = 'Elements of input list to segment_midpoints must be of class Triangle_intersection' |
---|
| 1233 | for segment in segments: |
---|
| 1234 | assert isinstance(segment, Triangle_intersection), msg |
---|
| 1235 | |
---|
[6174] | 1236 | midpoint = num.sum(num.array(segment.segment, num.Float))/2 |
---|
[5897] | 1237 | midpoints.append(midpoint) |
---|
| 1238 | |
---|
| 1239 | return midpoints |
---|
| 1240 | |
---|
| 1241 | |
---|
| 1242 | |
---|