Fix EMPTY result from GEOSOffsetCurve with distance 0 (#454)
[geos.git] / src / operation / buffer / BufferBuilder.cpp
blob77815d5ca987a045624b999d720af3f79b60c528
1 /**********************************************************************
2 * $Id$
4 * GEOS - Geometry Engine Open Source
5 * http://geos.refractions.net
7 * Copyright (C) 2009-2011 Sandro Santilli <strk@keybit.net>
8 * Copyright (C) 2008-2010 Safe Software Inc.
9 * Copyright (C) 2005-2007 Refractions Research Inc.
10 * Copyright (C) 2001-2002 Vivid Solutions Inc.
12 * This is free software; you can redistribute and/or modify it under
13 * the terms of the GNU Lesser General Public Licence as published
14 * by the Free Software Foundation.
15 * See the COPYING file for more information.
17 **********************************************************************
19 * Last port: operation/buffer/BufferBuilder.java r378 (JTS-1.12)
21 **********************************************************************/
23 #include <geos/geom/GeometryFactory.h>
24 #include <geos/geom/Location.h>
25 #include <geos/geom/Geometry.h>
26 #include <geos/geom/Polygon.h>
27 #include <geos/geom/GeometryCollection.h>
28 #include <geos/geom/LineString.h>
29 #include <geos/geom/MultiLineString.h>
30 #include <geos/operation/buffer/BufferBuilder.h>
31 #include <geos/operation/buffer/OffsetCurveBuilder.h>
32 #include <geos/operation/buffer/OffsetCurveSetBuilder.h>
33 #include <geos/operation/buffer/BufferSubgraph.h>
34 #include <geos/operation/buffer/SubgraphDepthLocater.h>
35 #include <geos/operation/overlay/OverlayOp.h>
36 #include <geos/operation/overlay/snap/SnapOverlayOp.h>
37 #include <geos/operation/overlay/PolygonBuilder.h>
38 #include <geos/operation/overlay/OverlayNodeFactory.h>
39 #include <geos/operation/linemerge/LineMerger.h>
40 #include <geos/algorithm/LineIntersector.h>
41 #include <geos/noding/IntersectionAdder.h>
42 #include <geos/noding/SegmentString.h>
43 #include <geos/noding/MCIndexNoder.h>
44 #include <geos/noding/NodedSegmentString.h>
45 #include <geos/geomgraph/Position.h>
46 #include <geos/geomgraph/PlanarGraph.h>
47 #include <geos/geomgraph/Label.h>
48 #include <geos/geomgraph/Node.h>
49 #include <geos/geomgraph/Edge.h>
50 #include <geos/util/GEOSException.h>
51 #include <geos/io/WKTWriter.h> // for debugging
52 #include <geos/util/IllegalArgumentException.h>
53 #include <geos/profiler.h>
55 #include <cassert>
56 #include <vector>
57 #include <iomanip>
58 #include <algorithm>
59 #include <iostream>
61 // Debug single sided buffer
62 //#define GEOS_DEBUG_SSB 1
64 #ifndef GEOS_DEBUG
65 #define GEOS_DEBUG 0
66 #endif
68 #ifndef JTS_DEBUG
69 #define JTS_DEBUG 0
70 #endif
72 //using namespace std;
73 using namespace geos::geom;
74 using namespace geos::geomgraph;
75 using namespace geos::noding;
76 using namespace geos::algorithm;
77 using namespace geos::operation::overlay;
78 using namespace geos::operation::linemerge;
80 namespace {
82 // Debug routine
83 template <class Iterator>
84 std::auto_ptr<Geometry>
85 convertSegStrings(const GeometryFactory* fact, Iterator it, Iterator et)
87 std::vector<Geometry*> lines;
88 while(it != et) {
89 const SegmentString* ss = *it;
90 LineString* line = fact->createLineString(ss->getCoordinates());
91 lines.push_back(line);
92 ++it;
94 return std::auto_ptr<Geometry>(fact->buildGeometry(lines));
99 namespace geos {
100 namespace operation { // geos.operation
101 namespace buffer { // geos.operation.buffer
103 #if PROFILE
104 static Profiler *profiler = Profiler::instance();
105 #endif
108 * Compute the change in depth as an edge is crossed from R to L
111 BufferBuilder::depthDelta(Label *label)
113 int lLoc=label->getLocation(0, Position::LEFT);
114 int rLoc=label->getLocation(0, Position::RIGHT);
115 if (lLoc== Location::INTERIOR && rLoc== Location::EXTERIOR)
116 return 1;
117 else if (lLoc== Location::EXTERIOR && rLoc== Location::INTERIOR)
118 return -1;
119 return 0;
122 //static CGAlgorithms rCGA;
123 //CGAlgorithms *BufferBuilder::cga=&rCGA;
125 BufferBuilder::~BufferBuilder()
127 delete li; // could be NULL
128 delete intersectionAdder;
129 //delete edgeList;
130 for (size_t i=0; i<newLabels.size(); i++)
131 delete newLabels[i];
134 /*public*/
135 Geometry*
136 BufferBuilder::bufferLineSingleSided( const Geometry* g, double distance,
137 bool leftSide )
139 // Returns the line used to create a single-sided buffer.
140 // Input requirement: Must be a LineString.
141 const LineString* l = dynamic_cast< const LineString* >( g );
142 if ( !l ) throw util::IllegalArgumentException("BufferBuilder::bufferLineSingleSided only accept linestrings");
144 // Nothing to do for a distance of zero
145 if ( distance == 0 ) return g->clone();
147 // Get geometry factory and precision model.
148 const PrecisionModel* precisionModel = workingPrecisionModel;
149 if ( !precisionModel ) precisionModel = l->getPrecisionModel();
151 assert( precisionModel );
152 assert( l );
154 geomFact = l->getFactory();
156 // First, generate the two-sided buffer using a butt-cap.
157 BufferParameters modParams = bufParams;
158 modParams.setEndCapStyle(BufferParameters::CAP_FLAT);
159 Geometry* buf = 0;
161 // This is a (temp?) hack to workaround the fact that
162 // BufferBuilder BufferParamaters are immutable after
163 // construction, while we want to force the end cap
164 // style to FLAT for single-sided buffering
166 BufferBuilder tmp(modParams);
167 buf = tmp.buffer( l, distance );
170 // Create MultiLineStrings from this polygon.
171 Geometry* bufLineString = buf->getBoundary();
173 #ifdef GEOS_DEBUG_SSB
174 std::cerr << "input|" << *l << std::endl;
175 std::cerr << "buffer|" << *bufLineString << std::endl;
176 #endif
178 // Then, get the raw (i.e. unnoded) single sided offset curve.
179 OffsetCurveBuilder curveBuilder( precisionModel, modParams );
180 std::vector< CoordinateSequence* > lineList;
182 std::auto_ptr< CoordinateSequence > coords ( g->getCoordinates() );
183 curveBuilder.getSingleSidedLineCurve( coords.get(), distance,
184 lineList, leftSide, !leftSide );
185 coords.reset();
187 // Create a SegmentString from these coordinates.
188 SegmentString::NonConstVect curveList;
189 for ( unsigned int i = 0; i < lineList.size(); ++i )
191 CoordinateSequence* seq = lineList[i];
193 SegmentString* ss = new NodedSegmentString(seq, NULL);
194 curveList.push_back( ss );
197 // Node these SegmentStrings.
198 Noder* noder = getNoder( precisionModel );
199 noder->computeNodes( &curveList );
200 SegmentString::NonConstVect* nodedEdges = noder->getNodedSubstrings();
202 // Create a geometry out of the noded substrings.
203 std::vector< Geometry* >* singleSidedNodedEdges =
204 new std::vector< Geometry * >();
205 for ( unsigned int i = 0, n = nodedEdges->size(); i < n; ++i )
207 SegmentString* ss = ( *nodedEdges )[i];
209 Geometry* tmp = geomFact->createLineString(
210 ss->getCoordinates()->clone()
212 singleSidedNodedEdges->push_back( tmp );
215 if ( nodedEdges != &curveList ) delete nodedEdges;
217 for (size_t i=0, n=curveList.size(); i<n; ++i) delete curveList[i];
218 curveList.clear();
220 for (size_t i=0, n=lineList.size(); i<n; ++i) delete lineList[i];
221 lineList.clear();
224 Geometry* singleSided = geomFact->createMultiLineString(
225 singleSidedNodedEdges );
227 #ifdef GEOS_DEBUG_SSB
228 std::cerr << "edges|" << *singleSided << std::endl;
229 #endif
231 // Use the boolean operation intersect to obtain the line segments lying
232 // on both the butt-cap buffer and this multi-line.
233 //Geometry* intersectedLines = singleSided->intersection( bufLineString );
234 // NOTE: we use Snapped overlay because the actual buffer boundary might
235 // diverge from original offset curves due to the addition of
236 // intersections with caps and joins curves
237 using geos::operation::overlay::snap::SnapOverlayOp;
238 Geometry* intersectedLines = SnapOverlayOp::overlayOp(*singleSided, *bufLineString, OverlayOp::opINTERSECTION).release();
240 #ifdef GEOS_DEBUG_SSB
241 std::cerr << "intersection" << "|" << *intersectedLines << std::endl;
242 #endif
244 // Merge result lines together.
245 LineMerger lineMerge;
246 lineMerge.add( intersectedLines );
247 std::auto_ptr< std::vector< LineString* > > mergedLines (
248 lineMerge.getMergedLineStrings() );
250 // Convert the result into a std::vector< Geometry* >.
251 std::vector< Geometry* >* mergedLinesGeom = new std::vector< Geometry* >();
252 const Coordinate& startPoint = l->getCoordinatesRO()->front();
253 const Coordinate& endPoint = l->getCoordinatesRO()->back();
254 while( !mergedLines->empty() )
256 // Remove end points if they are a part of the original line to be
257 // buffered.
258 CoordinateSequence::AutoPtr coords(mergedLines->back()->getCoordinates());
259 if ( NULL != coords.get() )
261 // Use 98% of the buffer width as the point-distance requirement - this
262 // is to ensure that the point that is "distance" +/- epsilon is not
263 // included.
265 // Let's try and estimate a more accurate bound instead of just assuming
266 // 98%. With 98%, the episilon grows as the buffer distance grows,
267 // so that at large distance, artifacts may skip through this filter
268 // Let the length of the line play a factor in the distance, which is still
269 // going to be bounded by 98%. Take 10% of the length of the line from the buffer distance
270 // to try and minimize any artifacts.
271 const double ptDistAllowance = (std::max)(distance - l->getLength()*0.1, distance * 0.98);
272 // Use 102% of the buffer width as the line-length requirement - this
273 // is to ensure that line segments that is length "distance" +/-
274 // epsilon is removed.
275 const double segLengthAllowance = 1.02 * distance;
277 // Clean up the front of the list.
278 // Loop until the line's end is not inside the buffer width from
279 // the startPoint.
280 while ( coords->size() > 1 &&
281 coords->front().distance( startPoint ) < ptDistAllowance )
283 // Record the end segment length.
284 double segLength = coords->front().distance( ( *coords )[1] );
285 // Stop looping if there are no more points, or if the segment
286 // length is larger than the buffer width.
287 if ( coords->size() <= 1 || segLength > segLengthAllowance )
289 break;
291 // If the first point is less than buffer width away from the
292 // reference point, then delete the point.
293 coords->deleteAt( 0 );
295 while ( coords->size() > 1 &&
296 coords->front().distance( endPoint ) < ptDistAllowance )
298 double segLength = coords->front().distance( ( *coords )[1] );
299 if ( coords->size() <= 1 || segLength > segLengthAllowance )
301 break;
303 coords->deleteAt( 0 );
306 // Clean up the back of the list.
307 while ( coords->size() > 1 &&
308 coords->back().distance( startPoint ) < ptDistAllowance )
310 double segLength = coords->back().distance(
311 ( *coords )[coords->size()-2] );
312 if ( coords->size() <= 1 || segLength > segLengthAllowance )
314 break;
316 coords->deleteAt( coords->size()-1 );
318 while ( coords->size() > 1 &&
319 coords->back().distance( endPoint ) < ptDistAllowance )
321 double segLength = coords->back().distance(
322 ( *coords )[coords->size()-2] );
323 if ( coords->size() <= 1 || segLength > segLengthAllowance )
325 break;
327 coords->deleteAt( coords->size()-1 );
330 // Add the coordinates to the resultant line string.
331 if ( coords->size() > 1 )
333 mergedLinesGeom->push_back( geomFact->createLineString( coords.release() ) );
337 geomFact->destroyGeometry( mergedLines->back() );
338 mergedLines->pop_back();
341 // Clean up.
342 if ( noder != workingNoder ) delete noder;
343 geomFact->destroyGeometry( buf );
344 geomFact->destroyGeometry( bufLineString );
345 geomFact->destroyGeometry( singleSided );
346 geomFact->destroyGeometry( intersectedLines );
348 if ( mergedLinesGeom->size() > 1 )
350 return geomFact->createMultiLineString( mergedLinesGeom );
352 else if ( mergedLinesGeom->size() == 1 )
355 Geometry* single = (*mergedLinesGeom)[0];
356 delete mergedLinesGeom;
357 return single;
359 else
361 delete mergedLinesGeom;
362 return geomFact->createLineString();
366 /*public*/
367 Geometry*
368 BufferBuilder::buffer(const Geometry *g, double distance)
369 // throw(GEOSException *)
371 const PrecisionModel *precisionModel=workingPrecisionModel;
372 if (precisionModel==NULL)
373 precisionModel=g->getPrecisionModel();
375 assert(precisionModel);
376 assert(g);
378 // factory must be the same as the one used by the input
379 geomFact=g->getFactory();
381 OffsetCurveBuilder curveBuilder(precisionModel, bufParams);
382 OffsetCurveSetBuilder curveSetBuilder(*g, distance, curveBuilder);
384 std::vector<SegmentString*>& bufferSegStrList=curveSetBuilder.getCurves();
386 #if GEOS_DEBUG
387 std::cerr << "OffsetCurveSetBuilder got " << bufferSegStrList.size()
388 << " curves" << std::endl;
389 #endif
390 // short-circuit test
391 if (bufferSegStrList.size()<=0) {
392 return createEmptyResultGeometry();
395 #if GEOS_DEBUG
396 std::cerr<<"BufferBuilder::buffer computing NodedEdges"<<std::endl;
397 #endif
399 computeNodedEdges(bufferSegStrList, precisionModel);
400 // NOTE: bufferSegStrList should not be needed anymore from now on
402 #if GEOS_DEBUG > 1
403 std::cerr << std::endl << edgeList << std::endl;
404 #endif
406 Geometry* resultGeom=NULL;
407 std::auto_ptr< std::vector<Geometry*> > resultPolyList;
408 std::vector<BufferSubgraph*> subgraphList;
410 try {
411 PlanarGraph graph(OverlayNodeFactory::instance());
412 graph.addEdges(edgeList.getEdges());
414 createSubgraphs(&graph, subgraphList);
415 #if GEOS_DEBUG
416 std::cerr<<"Created "<<subgraphList.size()<<" subgraphs"<<std::endl;
417 #if GEOS_DEBUG > 1
418 for (size_t i=0, n=subgraphList.size(); i<n; i++)
419 std::cerr << std::setprecision(10) << *(subgraphList[i]) << std::endl;
420 #endif
421 #endif
422 PolygonBuilder polyBuilder(geomFact);
423 buildSubgraphs(subgraphList, polyBuilder);
424 resultPolyList.reset( polyBuilder.getPolygons() );
425 #if GEOS_DEBUG
426 std::cerr << "PolygonBuilder got " << resultPolyList->size()
427 << " polygons" << std::endl;
428 #if GEOS_DEBUG > 1
429 for (size_t i=0, n=resultPolyList->size(); i<n; i++)
430 std::cerr << (*resultPolyList)[i]->toString() << std::endl;
431 #endif
432 #endif
433 // just in case ...
434 if ( resultPolyList->empty() )
436 for (size_t i=0, n=subgraphList.size(); i<n; i++)
437 delete subgraphList[i];
438 return createEmptyResultGeometry();
441 // resultPolyList ownership transferred here
442 resultGeom=geomFact->buildGeometry(resultPolyList.release());
443 } catch (const util::GEOSException& /* exc */) {
444 for (size_t i=0, n=subgraphList.size(); i<n; i++)
445 delete subgraphList[i];
446 throw;
449 for (size_t i=0, n=subgraphList.size(); i<n; i++)
450 delete subgraphList[i];
452 return resultGeom;
455 /*private*/
456 Noder*
457 BufferBuilder::getNoder(const PrecisionModel* pm)
459 // this doesn't change workingNoder precisionModel!
460 if (workingNoder != NULL) return workingNoder;
462 // otherwise use a fast (but non-robust) noder
464 if ( li ) // reuse existing IntersectionAdder and LineIntersector
466 li->setPrecisionModel(pm);
467 assert(intersectionAdder!=NULL);
469 else
471 li = new LineIntersector(pm);
472 intersectionAdder = new IntersectionAdder(*li);
475 MCIndexNoder* noder = new MCIndexNoder(intersectionAdder);
477 #if 0
478 /* CoordinateArraySequence.cpp:84:
479 * virtual const geos::Coordinate& geos::CoordinateArraySequence::getAt(size_t) const:
480 * Assertion `pos<vect->size()' failed.
482 //Noder* noder = new snapround::SimpleSnapRounder(*pm);
484 Noder* noder = new IteratedNoder(pm);
486 Noder noder = new SimpleSnapRounder(pm);
487 Noder noder = new MCIndexSnapRounder(pm);
488 Noder noder = new ScaledNoder(
489 new MCIndexSnapRounder(new PrecisionModel(1.0)),
490 pm.getScale());
491 #endif
493 return noder;
497 /* private */
498 void
499 BufferBuilder::computeNodedEdges(SegmentString::NonConstVect& bufferSegStrList,
500 const PrecisionModel *precisionModel) // throw(GEOSException)
502 Noder* noder = getNoder( precisionModel );
504 #if JTS_DEBUG
505 geos::io::WKTWriter wktWriter; wktWriter.setTrim(true);
506 std::cerr << "before noding: "
507 << wktWriter.write(
508 convertSegStrings(geomFact, bufferSegStrList.begin(),
509 bufferSegStrList.end()).get()
510 ) << std::endl;
511 #endif
513 noder->computeNodes(&bufferSegStrList);
515 SegmentString::NonConstVect* nodedSegStrings = \
516 noder->getNodedSubstrings();
518 #if JTS_DEBUG
519 std::cerr << "after noding: "
520 << wktWriter.write(
521 convertSegStrings(geomFact, bufferSegStrList.begin(),
522 bufferSegStrList.end()).get()
523 ) << std::endl;
524 #endif
527 for (SegmentString::NonConstVect::iterator
528 i=nodedSegStrings->begin(), e=nodedSegStrings->end();
529 i!=e;
530 ++i)
532 SegmentString* segStr = *i;
533 const Label* oldLabel = static_cast<const Label*>(segStr->getData());
535 CoordinateSequence* cs = CoordinateSequence::removeRepeatedPoints(segStr->getCoordinates());
536 if ( cs->size() < 2 )
538 delete cs; // we need to take care of the memory here as cs is a new sequence
539 return; // don't insert collapsed edges
541 // we need to clone SegmentString coordinates
542 // as Edge will take ownership of them
543 // TODO: find a way to transfer ownership instead
544 // Who will own the edge ? FIXME: find out and handle that!
545 Edge* edge = new Edge(cs, new Label(*oldLabel));
547 // will take care of the Edge ownership
548 insertUniqueEdge(edge);
551 if ( nodedSegStrings != &bufferSegStrList )
553 delete nodedSegStrings;
556 if ( noder != workingNoder ) delete noder;
559 /*private*/
560 void
561 BufferBuilder::insertUniqueEdge(Edge *e)
563 //<FIX> MD 8 Oct 03 speed up identical edge lookup
564 // fast lookup
565 Edge *existingEdge=edgeList.findEqualEdge(e);
566 // If an identical edge already exists, simply update its label
567 if (existingEdge != NULL) {
568 Label *existingLabel=existingEdge->getLabel();
569 Label *labelToMerge=e->getLabel();
571 // check if new edge is in reverse direction to existing edge
572 // if so, must flip the label before merging it
573 if (! existingEdge->isPointwiseEqual(e))
575 labelToMerge=new Label(*(e->getLabel()));
576 labelToMerge->flip();
577 newLabels.push_back(labelToMerge);
579 existingLabel->merge(*labelToMerge);
580 // compute new depth delta of sum of edges
581 int mergeDelta=depthDelta(labelToMerge);
582 int existingDelta=existingEdge->getDepthDelta();
583 int newDelta=existingDelta + mergeDelta;
584 existingEdge->setDepthDelta(newDelta);
586 // we have memory release responsibility
587 delete e;
589 } else { // no matching existing edge was found
591 // add this new edge to the list of edges in this graph
592 edgeList.add(e);
594 e->setDepthDelta(depthDelta(e->getLabel()));
598 bool BufferSubgraphGT(BufferSubgraph *first, BufferSubgraph *second) {
599 if (first->compareTo(second)>0)
600 return true;
601 else
602 return false;
605 /*private*/
606 void
607 BufferBuilder::createSubgraphs(PlanarGraph *graph, std::vector<BufferSubgraph*>& subgraphList)
609 std::vector<Node*> nodes;
610 graph->getNodes(nodes);
611 for (size_t i=0, n=nodes.size(); i<n; i++) {
612 Node *node=nodes[i];
613 if (!node->isVisited()) {
614 BufferSubgraph *subgraph=new BufferSubgraph();
615 subgraph->create(node);
616 subgraphList.push_back(subgraph);
621 * Sort the subgraphs in descending order of their rightmost coordinate
622 * This ensures that when the Polygons for the subgraphs are built,
623 * subgraphs for shells will have been built before the subgraphs for
624 * any holes they contain
626 std::sort(subgraphList.begin(), subgraphList.end(), BufferSubgraphGT);
629 /*private*/
630 void
631 BufferBuilder::buildSubgraphs(const std::vector<BufferSubgraph*>& subgraphList,
632 PolygonBuilder& polyBuilder)
635 #if GEOS_DEBUG
636 std::cerr << __FUNCTION__ << " got " << subgraphList.size() << " subgraphs" << std::endl;
637 #endif
638 std::vector<BufferSubgraph*> processedGraphs;
639 for (size_t i=0, n=subgraphList.size(); i<n; i++)
641 BufferSubgraph *subgraph=subgraphList[i];
642 Coordinate *p=subgraph->getRightmostCoordinate();
643 assert(p);
645 #if GEOS_DEBUG
646 std::cerr << " " << i << ") Subgraph[" << subgraph << "]" << std::endl;
647 std::cerr << " rightmost Coordinate " << *p;
648 #endif
649 SubgraphDepthLocater locater(&processedGraphs);
650 #if GEOS_DEBUG
651 std::cerr << " after SubgraphDepthLocater processedGraphs contain "
652 << processedGraphs.size()
653 << " elements" << std::endl;
654 #endif
655 int outsideDepth=locater.getDepth(*p);
656 #if GEOS_DEBUG
657 std::cerr << " Depth of rightmost coordinate: " << outsideDepth << std::endl;
658 #endif
659 subgraph->computeDepth(outsideDepth);
660 subgraph->findResultEdges();
661 #if GEOS_DEBUG
662 std::cerr << " after computeDepth and findResultEdges subgraph contain:" << std::endl
663 << " " << subgraph->getDirectedEdges()->size() << " DirecteEdges " << std::endl
664 << " " << subgraph->getNodes()->size() << " Nodes " << std::endl;
665 #endif
666 processedGraphs.push_back(subgraph);
667 #if GEOS_DEBUG
668 std::cerr << " added " << subgraph << " to processedGraphs, new size is "
669 << processedGraphs.size() << std::endl;
670 #endif
671 polyBuilder.add(subgraph->getDirectedEdges(), subgraph->getNodes());
675 /*private*/
676 geom::Geometry*
677 BufferBuilder::createEmptyResultGeometry() const
679 geom::Geometry* emptyGeom = geomFact->createPolygon(NULL, NULL);
680 return emptyGeom;
683 } // namespace geos.operation.buffer
684 } // namespace geos.operation
685 } // namespace geos
687 /**********************************************************************
688 * $Log$
689 * Revision 1.56 2006/06/12 11:29:23 strk
690 * unsigned int => size_t
692 * Revision 1.55 2006/05/04 09:16:58 strk
693 * Added JTS debugging, for comparison with JTS
695 * Revision 1.54 2006/03/24 09:25:02 strk
696 * Bugs #77 and #76: missing <algorithm>
698 * Revision 1.53 2006/03/22 18:12:32 strk
699 * indexChain.h header split.
701 * Revision 1.52 2006/03/17 13:24:59 strk
702 * opOverlay.h header splitted. Reduced header inclusions in operation/overlay implementation files. ElevationMatrixFilter code moved from own file to ElevationMatrix.cpp (ideally a class-private).
704 * Revision 1.51 2006/03/15 18:57:10 strk
705 * cleanups in DEBUG lines
707 * Revision 1.50 2006/03/15 13:03:01 strk
708 * removed leftover debugging line
710 * Revision 1.49 2006/03/15 11:42:54 strk
711 * more debugging lines, with two levels of debugging handled
713 * Revision 1.48 2006/03/14 16:08:21 strk
714 * changed buildSubgraphs signature to use refs rather then pointers, made it const-correct. Reduced heap allocations in createSubgraphs()
716 * Revision 1.47 2006/03/14 14:16:52 strk
717 * operator<< for BufferSubgraph, more debugging calls
719 * Revision 1.46 2006/03/14 12:55:56 strk
720 * Headers split: geomgraphindex.h, nodingSnapround.h
722 * Revision 1.45 2006/03/14 00:19:40 strk
723 * opBuffer.h split, streamlined headers in some (not all) files in operation/buffer/
725 * Revision 1.44 2006/03/10 10:44:53 strk
726 * Unreferenced exception objects cleanup (#52)
728 * Revision 1.43 2006/03/06 19:40:47 strk
729 * geos::util namespace. New GeometryCollection::iterator interface, many cleanups.
731 * Revision 1.42 2006/03/03 10:46:21 strk
732 * Removed 'using namespace' from headers, added missing headers in .cpp files, removed useless includes in headers (bug#46)
734 * Revision 1.41 2006/03/02 12:12:01 strk
735 * Renamed DEBUG macros to GEOS_DEBUG, all wrapped in #ifndef block to allow global override (bug#43)
737 * Revision 1.40 2006/02/28 17:44:27 strk
738 * Added a check in SegmentNode::addSplitEdge to prevent attempts
739 * to build SegmentString with less then 2 points.
740 * This is a temporary fix for the buffer.xml assertion failure, temporary
741 * as Martin Davis review would really be needed there.
743 * Revision 1.39 2006/02/28 14:34:05 strk
744 * Added many assertions and debugging output hunting for a bug in BufferOp
746 * Revision 1.38 2006/02/23 20:05:21 strk
747 * Fixed bug in MCIndexNoder constructor making memory checker go crazy, more
748 * doxygen-friendly comments, miscellaneous cleanups
750 * Revision 1.37 2006/02/23 11:54:20 strk
751 * - MCIndexPointSnapper
752 * - MCIndexSnapRounder
753 * - SnapRounding BufferOp
754 * - ScaledNoder
755 * - GEOSException hierarchy cleanups
756 * - SpatialIndex memory-friendly query interface
757 * - GeometryGraph::getBoundaryNodes memory-friendly
758 * - NodeMap::getBoundaryNodes memory-friendly
759 * - Cleanups in geomgraph::Edge
760 * - Added an XML test for snaprounding buffer (shows leaks, working on it)
762 * Revision 1.36 2006/02/19 19:46:49 strk
763 * Packages <-> namespaces mapping for most GEOS internal code (uncomplete, but working). Dir-level libs for index/ subdirs.
765 * Revision 1.35 2006/02/18 21:08:09 strk
766 * - new CoordinateSequence::applyCoordinateFilter method (slow but useful)
767 * - SegmentString::getCoordinates() doesn't return a clone anymore.
768 * - SegmentString::getCoordinatesRO() obsoleted.
769 * - SegmentString constructor does not promises constness of passed
770 * CoordinateSequence anymore.
771 * - NEW ScaledNoder class
772 * - Stubs for MCIndexPointSnapper and MCIndexSnapRounder
773 * - Simplified internal interaces of OffsetCurveBuilder and OffsetCurveSetBuilder
775 * Revision 1.34 2006/02/14 13:28:26 strk
776 * New SnapRounding code ported from JTS-1.7 (not complete yet).
777 * Buffer op optimized by using new snaprounding code.
778 * Leaks fixed in XMLTester.
780 * Revision 1.33 2006/02/09 15:52:47 strk
781 * GEOSException derived from std::exception; always thrown and cought by const ref.
783 * Revision 1.32 2006/02/08 17:18:28 strk
784 * - New WKTWriter::toLineString and ::toPoint convenience methods
785 * - New IsValidOp::setSelfTouchingRingFormingHoleValid method
786 * - New Envelope::centre()
787 * - New Envelope::intersection(Envelope)
788 * - New Envelope::expandBy(distance, [ydistance])
789 * - New LineString::reverse()
790 * - New MultiLineString::reverse()
791 * - New Geometry::buffer(distance, quadSeg, endCapStyle)
792 * - Obsoleted toInternalGeometry/fromInternalGeometry
793 * - More const-correctness in Buffer "package"
795 * Revision 1.31 2006/01/31 19:07:34 strk
796 * - Renamed DefaultCoordinateSequence to CoordinateArraySequence.
797 * - Moved GetNumGeometries() and GetGeometryN() interfaces
798 * from GeometryCollection to Geometry class.
799 * - Added getAt(int pos, Coordinate &to) funtion to CoordinateSequence class.
800 * - Reworked automake scripts to produce a static lib for each subdir and
801 * then link all subsystem's libs togheter
802 * - Moved C-API in it's own top-level dir capi/
803 * - Moved source/bigtest and source/test to tests/bigtest and test/xmltester
804 * - Fixed PointLocator handling of LinearRings
805 * - Changed CoordinateArrayFilter to reduce memory copies
806 * - Changed UniqueCoordinateArrayFilter to reduce memory copies
807 * - Added CGAlgorithms::isPointInRing() version working with
808 * Coordinate::ConstVect type (faster!)
809 * - Ported JTS-1.7 version of ConvexHull with big attention to
810 * memory usage optimizations.
811 * - Improved XMLTester output and user interface
812 * - geos::geom::util namespace used for geom/util stuff
813 * - Improved memory use in geos::geom::util::PolygonExtractor
814 * - New ShortCircuitedGeometryVisitor class
815 * - New operation/predicate package
817 * Revision 1.30 2005/11/21 16:03:20 strk
819 * Coordinate interface change:
820 * Removed setCoordinate call, use assignment operator
821 * instead. Provided a compile-time switch to
822 * make copy ctor and assignment operators non-inline
823 * to allow for more accurate profiling.
825 * Coordinate copies removal:
826 * NodeFactory::createNode() takes now a Coordinate reference
827 * rather then real value. This brings coordinate copies
828 * in the testLeaksBig.xml test from 654818 to 645991
829 * (tested in 2.1 branch). In the head branch Coordinate
830 * copies are 222198.
831 * Removed useless coordinate copies in ConvexHull
832 * operations
834 * STL containers heap allocations reduction:
835 * Converted many containers element from
836 * pointers to real objects.
837 * Made some use of .reserve() or size
838 * initialization when final container size is known
839 * in advance.
841 * Stateless classes allocations reduction:
842 * Provided ::instance() function for
843 * NodeFactories, to avoid allocating
844 * more then one (they are all
845 * stateless).
847 * HCoordinate improvements:
848 * Changed HCoordinate constructor by HCoordinates
849 * take reference rather then real objects.
850 * Changed HCoordinate::intersection to avoid
851 * a new allocation but rather return into a provided
852 * storage. LineIntersector changed to reflect
853 * the above change.
855 * Revision 1.29 2005/11/14 18:14:04 strk
856 * Reduced heap allocations made by TopologyLocation and Label objects.
857 * Enforced const-correctness on GraphComponent.
858 * Cleanups.
860 * Revision 1.28 2005/08/22 13:31:17 strk
861 * Fixed comparator functions used with STL sort() algorithm to
862 * implement StrictWeakOrdering semantic.
864 * Revision 1.27 2005/05/19 10:29:28 strk
865 * Removed some CGAlgorithms instances substituting them with direct calls
866 * to the static functions. Interfaces accepting CGAlgorithms pointers kept
867 * for backward compatibility but modified to make the argument optional.
868 * Fixed a small memory leak in OffsetCurveBuilder::getRingCurve.
869 * Inlined some smaller functions encountered during bug hunting.
870 * Updated Copyright notices in the touched files.
872 * Revision 1.26 2005/02/01 13:44:59 strk
873 * More profiling labels.
875 * Revision 1.25 2004/12/08 13:54:43 strk
876 * gcc warnings checked and fixed, general cleanups.
878 * Revision 1.24 2004/11/04 19:08:07 strk
879 * Cleanups, initializers list, profiling.
881 **********************************************************************/