Inkscape
Vector Graphics Editor
Loading...
Searching...
No Matches
cola.h
Go to the documentation of this file.
1/*
2 * vim: ts=4 sw=4 et tw=0 wm=0
3 *
4 * libcola - A library providing force-directed network layout using the
5 * stress-majorization method subject to separation constraints.
6 *
7 * Copyright (C) 2006-2015 Monash University
8 *
9 * This library is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public
11 * License as published by the Free Software Foundation; either
12 * version 2.1 of the License, or (at your option) any later version.
13 * See the file LICENSE.LGPL distributed with the library.
14 *
15 * This library is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
18 *
19 * Author(s): Tim Dwyer
20 * Michael Wybrow
21*/
22
23#ifndef COLA_H
24#define COLA_H
25
26#include <utility>
27#include <iterator>
28#include <vector>
29#include <valarray>
30#include <algorithm>
31#include <cmath>
32#include <iostream>
33
35#include "libcola/cluster.h"
37#include "libcola/exceptions.h"
39
40namespace vpsc { class Rectangle; }
41namespace topology {
42 class ColaTopologyAddon;
43}
44namespace dialect {
45 class Graph;
46}
47
48
56namespace cola {
57
58class NonOverlapConstraints;
59class NonOverlapConstraintExemptions;
60
62typedef std::vector<unsigned> NodeIndexes;
63
65typedef std::vector<NodeIndexes> ListOfNodeIndexes;
66
68typedef std::pair<unsigned, unsigned> Edge;
69
72typedef std::vector<double> EdgeLengths;
73#define StandardEdgeLengths EdgeLengths()
74
78class Lock {
79public:
80 Lock() {}
88 Lock(unsigned id, double X, double Y) : id(id), x(X), y(Y) {
89 }
90 unsigned getID() const {
91 return id;
92 }
93 double pos(vpsc::Dim dim) const {
94 return dim==vpsc::HORIZONTAL?x:y;
95 }
96private:
97 unsigned id;
98 double x;
99 double y;
100};
102typedef std::vector<cola::Lock> Locks;
103
107class Resize {
108public:
122 Resize(unsigned id, double x, double y, double w, double h)
123 : id(id), target(x,x+w,y,y+h) {}
124 unsigned getID() const {
125 return id;
126 }
127 const vpsc::Rectangle* getTarget() const {
128 return &target;
129 }
130private:
131 unsigned id;
133};
135typedef std::vector<cola::Resize> Resizes;
136
137/*
138 * Setting a desired position for a node adds a term to the goal function
139 * drawing the node towards that desired position
140 */
142 unsigned id;
143 double x;
144 double y;
145 double weight;
146};
147typedef std::vector<cola::DesiredPosition> DesiredPositions;
148
149/*
150 * desired positions which should override those computed by applying forces
151 * are passed in for a set of nodes. The first entry is the Node->id, the
152 * second is the desired position.
153 */
154typedef std::pair<unsigned,double> DesiredPositionInDim;
155typedef std::vector<DesiredPositionInDim> DesiredPositionsInDim;
156
169public:
189
190// To prevent C++ objects from being destroyed in garbage collected languages
191// when the libraries are called from SWIG, we hide the declarations of the
192// destructors and prevent generation of default destructors.
193#ifndef SWIG
194 virtual ~PreIteration() {}
195#endif
196 virtual bool operator()() { return true; }
200private:
203};
204
217public:
219 TestConvergence(const double tol = 1e-4, const unsigned maxiterations = 100)
220 : tolerance(tol),
222 { reset(); }
223 virtual ~TestConvergence() {}
224
225public:
226 virtual bool operator()(const double new_stress, std::valarray<double> & X, std::valarray<double> & Y)
227 {
228 COLA_UNUSED(X);
229 COLA_UNUSED(Y);
230
231 iterations++;
232 //std::cout<<"iteration="<<iterations<<", old_stress="<<old_stress<<", new_stress="<<new_stress<<std::endl;
233 if (old_stress == DBL_MAX) {
234 old_stress = new_stress;
235 return iterations >= maxiterations;
236 }
237 // converged if relative decrease in stress falls below threshold
238 // or if stress increases (shouldn't happen for straight majorization)
239 bool converged =
240 (old_stress - new_stress) / (new_stress + 1e-10) < tolerance
242 old_stress = new_stress;
243 return converged;
244 }
245 void reset() {
246 old_stress = DBL_MAX;
247 iterations = 0;
248 }
249 const double tolerance;
250 const unsigned maxiterations;
251 unsigned iterations;
252};
253
271public:
294 std::vector<Edge> const & es,
296 const double idealLength,
297 EdgeLengths eLengths = StandardEdgeLengths,
298 TestConvergence *doneTest = nullptr,
299 PreIteration* preIteration=nullptr,
300 bool useNeighbourStress = false);
307 constrainedLayout = true;
308 this->ccs=ccs;
309 }
310
312 constrainedLayout = true;
314 for (size_t i = 0; i < ccs.size(); ++i) {
315 ccsp->push_back(ccs.at(i));
316 }
317 this->ccs=ccsp;
318 }
319
345 void setStickyNodes(const double stickyWeight,
346 std::valarray<double> const & startX,
347 std::valarray<double> const & startY);
348
352 void setScaling(bool scaling) {
353 this->scaling=scaling;
354 }
360 this->externalSolver=externalSolver;
361 }
368 void setAvoidOverlaps(bool horizontal = false) {
369 constrainedLayout = true;
370 this->avoidOverlaps = horizontal ? Horizontal : Both;
371 }
385 void setStraightenEdges(std::vector<straightener::Edge*>* straightenEdges,
386 double bendWeight = 0.01, double potBendWeight = 0.1,
387 bool xSkipping = true) {
388 for(std::vector<straightener::Edge*>::const_iterator e=straightenEdges->begin();
389 e!=straightenEdges->end();e++) {
390 (*e)->rerouteAround(boundingBoxes);
391 }
392 constrainedLayout = true;
393 this->xSkipping = xSkipping;
394 this->straightenEdges = straightenEdges;
395 this->bendWeight = bendWeight;
397 }
402 for(unsigned i=0;i<n;i++) {
403 boundingBoxes[i]->moveCentre(X[i],Y[i]);
404 }
405 }
406
409 {
410 delete done;
411 }
412
414 delete gpX;
415 delete gpY;
416 }
417 }
427 void run(bool x=true, bool y=true);
439 void runOnce(bool x=true, bool y=true);
440 void straighten(std::vector<straightener::Edge*>&, vpsc::Dim);
443 }
444 double computeStress();
445private:
446 double euclidean_distance(unsigned i, unsigned j) {
447 return sqrt(
448 (X[i] - X[j]) * (X[i] - X[j]) +
449 (Y[i] - Y[j]) * (Y[i] - Y[j]));
450 }
451 double compute_stress(std::valarray<double> const & Dij);
452 void majorize(std::valarray<double> const & Dij,GradientProjection* gp, std::valarray<double>& coords, std::valarray<double> const & startCoords);
453 void newton(std::valarray<double> const & Dij,GradientProjection* gp, std::valarray<double>& coords, std::valarray<double> const & startCoords);
454 unsigned n; //< number of nodes
455 //std::valarray<double> degrees;
456 std::valarray<double> lap2; //< graph laplacian
457 std::valarray<double> Q; //< quadratic terms matrix used in computations
458 std::valarray<double> Dij; //< all pairs shortest path distances
459 double tol; //< convergence tolerance
460 TestConvergence *done; //< functor used to determine if layout is finished
461 bool using_default_done; // Whether we allocated a default TestConvergence object.
462 PreIteration* preIteration; //< client can use this to create locks on nodes
463 vpsc::Rectangles boundingBoxes; //< node bounding boxes
464 /*
465 * stickyNodes controls whether nodes are attracted to their starting
466 * positions (at time of ConstrainedMajorizationLayout instantiation)
467 * stored in startX, startY
468 */
469 std::valarray<double> X, Y;
472 std::valarray<double> startX;
473 std::valarray<double> startY;
477 /*
478 * A cluster is a set of nodes that are somehow semantically grouped
479 * and should therefore be kept together a bit more tightly than, and
480 * preferably without overlapping, the rest of the graph.
481 *
482 * We achieve this by augmenting the L matrix with stronger attractive
483 * forces between all members of a cluster (other than the root)
484 * and by maintaining a (preferably convex) hull around those
485 * constituents which, using constraints and dummy variables, is
486 * prevented from overlapping other parts of the graph.
487 *
488 * Clusters are defined over the graph in a hierarchy starting with
489 * a single root cluster.
490 *
491 * Need to:
492 * - augment Lap matrix with intra cluster forces
493 * - compute convex hull of each cluster
494 * - from convex hull generate "StraightenEdges"
495 */
501 std::vector<straightener::Edge*>* straightenEdges;
502
504 /*
505 * determines whether we should leave some overlaps to be resolved
506 * vertically when generating straightening constraints in the x-dim
507 */
509 /*
510 * when using the gradient projection optimisation method, the following
511 * controls whether the problem should be preconditioned by affine scaling
512 */
514 /*
515 * if the Mosek quadratic programming environment is available it may be
516 * used to solve each iteration of stress majorization... slow but useful
517 * for testing
518 */
521};
522
524
526
532{
533 public:
537
539 {
540 }
541
542 virtual TopologyAddonInterface *clone(void) const
543 {
544 return new TopologyAddonInterface(*this);
545 }
546
547 virtual void freeAssociatedObjects(void)
548 {
549 }
550
551 virtual void handleResizes(const Resizes& resizeList, unsigned n,
552 std::valarray<double>& X, std::valarray<double>& Y,
554 vpsc::Rectangles& boundingBoxes,
555 cola::RootCluster* clusterHierarchy)
556 {
557 COLA_UNUSED(resizeList);
558 COLA_UNUSED(n);
559 COLA_UNUSED(X);
560 COLA_UNUSED(Y);
561 COLA_UNUSED(ccs);
562 COLA_UNUSED(boundingBoxes);
563 COLA_UNUSED(clusterHierarchy);
564 }
565 virtual void computePathLengths(unsigned short** G)
566 {
567 COLA_UNUSED(G);
568 }
569 virtual double computeStress(void) const
570 {
571 return 0;
572 }
573 virtual bool useTopologySolver(void) const
574 {
575 return false;
576 }
577 virtual void makeFeasible(bool generateNonOverlapConstraints,
578 vpsc::Rectangles& boundingBoxes,
579 cola::RootCluster* clusterHierarchy)
580 {
581 COLA_UNUSED(generateNonOverlapConstraints);
582 COLA_UNUSED(boundingBoxes);
583 COLA_UNUSED(clusterHierarchy);
584 }
585 virtual void moveTo(const vpsc::Dim dim, vpsc::Variables& vs,
586 vpsc::Constraints& cs, std::valarray<double> &coords,
587 cola::RootCluster* clusterHierarchy)
588 {
589 COLA_UNUSED(dim);
590 COLA_UNUSED(vs);
591 COLA_UNUSED(cs);
592 COLA_UNUSED(coords);
593 COLA_UNUSED(clusterHierarchy);
594 }
596 const vpsc::Dim dim, std::valarray<double>& g,
598 std::valarray<double> &coords,
599 DesiredPositionsInDim& des, double oldStress)
600 {
601 COLA_UNUSED(layout);
602 COLA_UNUSED(dim);
603 COLA_UNUSED(g);
604 COLA_UNUSED(vs);
605 COLA_UNUSED(cs);
606 COLA_UNUSED(coords);
607 COLA_UNUSED(des);
608 COLA_UNUSED(oldStress);
609 return 0;
610 }
611};
612
613
623public:
653 const vpsc::Rectangles& rs,
654 const std::vector<cola::Edge>& es,
655 const double idealLength,
656 const EdgeLengths& eLengths = StandardEdgeLengths,
657 TestConvergence* doneTest = nullptr,
658 PreIteration* preIteration = nullptr);
660
670 void run(bool x=true, bool y=true);
671
683 void runOnce(bool x=true, bool y=true);
684
691
707 void setAvoidNodeOverlaps(bool avoidOverlaps,
708 ListOfNodeIndexes listOfNodeGroups =
710
724
726
734 {
735 clusterHierarchy = hierarchy;
736 }
753 UnsatisfiableConstraintInfos *unsatisfiableX,
754 UnsatisfiableConstraintInfos *unsatisfiableY)
755 {
756 unsatisfiable.resize(2);
757 unsatisfiable[0]=unsatisfiableX;
758 unsatisfiable[1]=unsatisfiableY;
759 }
780 void makeFeasible(double xBorder=1, double yBorder=1);
781
795 void freeAssociatedObjects(void);
796
807 void outputInstanceToSVG(std::string filename = std::string());
808
820 void setUseNeighbourStress(bool useNeighbourStress);
821
836 std::vector<double> readLinearD(void);
837
857 std::vector<unsigned> readLinearG(void);
858
859 double computeStress() const;
860
861private:
862 unsigned n; // number of nodes
863 std::valarray<double> X, Y;
865 double applyForcesAndConstraints(const vpsc::Dim dim,const double oldStress);
866 double computeStepSize(const SparseMatrix& H, const std::valarray<double>& g,
867 const std::valarray<double>& d) const;
868 void computeDescentVectorOnBothAxes(const bool xaxis, const bool yaxis,
869 double stress, std::valarray<double>& x0, std::valarray<double>& x1);
870 void moveTo(const vpsc::Dim dim, std::valarray<double>& target);
871 double applyDescentVector(
872 const std::valarray<double>& d,
873 const std::valarray<double>& oldCoords,
874 std::valarray<double> &coords,
875 const double oldStress,
876 double stepsize
877 /*,topology::TopologyConstraints *s=nullptr*/);
879 const std::vector<Edge>& es, std::valarray<double> eLengths);
881 vpsc::Variables (&vs)[2]);
882 void handleResizes(const Resizes&);
883 void setPosition(std::valarray<double>& pos);
884 void moveBoundingBoxes();
885 bool noForces(double, double, unsigned) const;
886 void computeForces(const vpsc::Dim dim, SparseMap &H,
887 std::valarray<double> &g);
889 vpsc::Variables (&vars)[2], unsigned int& priority,
890 cola::NonOverlapConstraints *noc, Cluster *cluster,
891 cola::CompoundConstraints& idleConstraints);
892 std::vector<double> offsetDir(double minD);
893
894 void computeNeighbours(std::vector<Edge> es);
895 std::vector<std::vector<unsigned> > neighbours;
896 std::vector<std::vector<double> > neighbourLengths;
898 bool using_default_done; // Whether we allocated a default TestConvergence object.
901 double** D;
902 unsigned short** G;
903 double minD;
905
907 std::vector<UnsatisfiableConstraintInfos*> unsatisfiable;
911
917 const std::valarray<double> m_edge_lengths;
918
920
922 friend class dialect::Graph;
923};
924
927 std::string unsatinfo;
928};
929
948 bool preventOverlaps, int accept=0, unsigned debugLevel=0);
949
966 unsigned debugLevel=0);
967
968
971 std::vector<Edge> const & es,
972 RootCluster* clusterHierarchy,
973 const double idealLength,
974 bool useNeighbourStress = false
975 );
976
977/*
978 * find shortest path lengths from node s to all other nodes.
979 * @param s starting node
980 * @param n total number of nodes
981 * @param d n vector of path lengths
982 * @param es edge pairs
983 * @param eLengths edge weights
984 */
985void dijkstra(const unsigned s, const unsigned n, double* d,
986 const std::vector<cola::Edge>& es,
987 const std::valarray<double> & eLengths);
988
989#if 0
990void removeClusterOverlapFast(RootCluster& clusterHierarchy, vpsc::Rectangles& rs, Locks& locks);
991#endif
992
993void setupVarsAndConstraints(unsigned n, const CompoundConstraints& ccs,
994 const vpsc::Dim dim, vpsc::Rectangles& boundingBoxes,
995 RootCluster *clusterHierarchy,
997 std::valarray<double> &coords);
999 const DesiredPositionsInDim& des, std::valarray<double>& coords);
1000
1001}
1002#endif // COLA_H
AAF yaxis(AAF x, AAF y)
Definition aa.cpp:281
AAF xaxis(AAF x, AAF y)
Definition aa.cpp:269
Geom::IntRect bounds
Definition canvas.cpp:182
A cluster defines a hierarchical partitioning over the nodes which should be kept disjoint by the lay...
Definition cluster.h:51
Implements a constrained force-directed layout algorithm.
Definition cola.h:622
RootCluster * clusterHierarchy
Definition cola.h:912
void run(bool x=true, bool y=true)
Implements the main layout loop, taking descent steps until stress is no-longer significantly reduced...
Definition colafd.cpp:316
void computeDescentVectorOnBothAxes(const bool xaxis, const bool yaxis, double stress, std::valarray< double > &x0, std::valarray< double > &x1)
Definition colafd.cpp:297
void freeAssociatedObjects(void)
A convenience method that can be called from Java to free the memory of nodes (Rectangles),...
Definition colafd.cpp:916
DesiredPositions * desiredPositions
Definition cola.h:909
bool m_generateNonOverlapConstraints
Definition cola.h:915
void computeNeighbours(std::vector< Edge > es)
Definition colafd.cpp:171
void setConstraints(const cola::CompoundConstraints &ccs)
Specify a set of compound constraints to apply to the layout.
Definition colafd.cpp:189
std::vector< std::vector< unsigned > > neighbours
Definition cola.h:895
void computePathLengths(const std::vector< Edge > &es, std::valarray< double > eLengths)
Definition colafd.cpp:227
TestConvergence * done
Definition cola.h:897
std::vector< double > offsetDir(double minD)
Definition colafd.cpp:1199
double applyDescentVector(const std::valarray< double > &d, const std::valarray< double > &oldCoords, std::valarray< double > &coords, const double oldStress, double stepsize)
Definition colafd.cpp:1172
void setPosition(std::valarray< double > &pos)
Definition colafd.cpp:286
TopologyAddonInterface * topologyAddon
Definition cola.h:906
std::vector< std::vector< double > > neighbourLengths
Definition cola.h:896
double applyForcesAndConstraints(const vpsc::Dim dim, const double oldStress)
Definition colafd.cpp:1105
void setUnsatisfiableConstraintInfo(UnsatisfiableConstraintInfos *unsatisfiableX, UnsatisfiableConstraintInfos *unsatisfiableY)
Register to receive information about unsatisfiable constraints.
Definition cola.h:752
vpsc::Rectangles boundingBoxes
Definition cola.h:864
TopologyAddonInterface * getTopology(void)
Definition colafd.cpp:951
std::vector< unsigned > readLinearG(void)
Retrieve a copy of the "G matrix" computed by the computePathLengths method, linearised as a vector.
Definition colafd.cpp:159
double computeStepSize(const SparseMatrix &H, const std::valarray< double > &g, const std::valarray< double > &d) const
Definition colafd.cpp:1292
void setTopology(TopologyAddonInterface *topology)
Set an addon for doing topology preserving layout.
Definition colafd.cpp:944
NonOverlapConstraintExemptions * m_nonoverlap_exemptions
Definition cola.h:919
void setDesiredPositions(DesiredPositions *desiredPositions)
Definition colafd.cpp:206
std::valarray< double > Y
Definition cola.h:863
double computeStress() const
Definition colafd.cpp:1314
std::valarray< double > X
Definition cola.h:863
void computeForces(const vpsc::Dim dim, SparseMap &H, std::valarray< double > &g)
Definition colafd.cpp:1225
void handleResizes(const Resizes &)
Definition colafd.cpp:1052
std::vector< double > readLinearD(void)
Retrieve a copy of the "D matrix" computed by the computePathLengths method, linearised as a vector.
Definition colafd.cpp:147
std::vector< UnsatisfiableConstraintInfos * > unsatisfiable
Definition cola.h:907
const std::valarray< double > m_edge_lengths
Definition cola.h:917
bool noForces(double, double, unsigned) const
void moveTo(const vpsc::Dim dim, std::valarray< double > &target)
Definition colafd.cpp:1063
PseudoRandom random
Definition cola.h:904
void makeFeasible(double xBorder=1, double yBorder=1)
Finds a feasible starting position for nodes that satisfies the given constraints.
Definition colafd.cpp:604
void setAvoidNodeOverlaps(bool avoidOverlaps, ListOfNodeIndexes listOfNodeGroups=ListOfNodeIndexes())
Specifies whether non-overlap constraints should be automatically generated between all nodes,...
Definition colafd.cpp:194
cola::CompoundConstraints ccs
Definition cola.h:900
void generateNonOverlapAndClusterCompoundConstraints(vpsc::Variables(&vs)[2])
Definition colafd.cpp:535
void recGenerateClusterVariablesAndConstraints(vpsc::Variables(&vars)[2], unsigned int &priority, cola::NonOverlapConstraints *noc, Cluster *cluster, cola::CompoundConstraints &idleConstraints)
Definition colafd.cpp:416
friend class topology::ColaTopologyAddon
Definition cola.h:921
void outputInstanceToSVG(std::string filename=std::string())
Generates an SVG file containing debug output and code that can be used to regenerate the instance.
Definition colafd.cpp:1369
cola::CompoundConstraints extraConstraints
Definition cola.h:910
PreIteration * preIteration
Definition cola.h:899
unsigned short ** G
Definition cola.h:902
void runOnce(bool x=true, bool y=true)
Same as run(), but only applies one iteration.
Definition colafd.cpp:386
void setUseNeighbourStress(bool useNeighbourStress)
Specifies whether neighbour stress should be used.
Definition colafd.cpp:201
friend class dialect::Graph
Definition cola.h:922
void setClusterHierarchy(RootCluster *hierarchy)
Specifies an optional hierarchy for clustering nodes.
Definition cola.h:733
Implements the Constrained Majorization graph layout algorithm (deprecated).
Definition cola.h:270
std::valarray< double > Q
Definition cola.h:457
GradientProjection * gpY
Definition cola.h:497
void setAvoidOverlaps(bool horizontal=false)
At each iteration of layout, generate constraints to avoid overlaps.
Definition cola.h:368
std::valarray< double > X
Definition cola.h:469
void runOnce(bool x=true, bool y=true)
Same as run(), but only applies one iteration.
Definition cola.cpp:397
std::valarray< double > lap2
Definition cola.h:456
std::valarray< double > Dij
Definition cola.h:458
void setConstraintsVector(cola::CompoundConstraints &ccs)
Definition cola.h:311
void setExternalSolver(bool externalSolver)
Says that the Mosek optimisation library should be used to solve the quadratic programs rather than t...
Definition cola.h:359
double euclidean_distance(unsigned i, unsigned j)
Definition cola.h:446
cola::CompoundConstraints * ccs
Definition cola.h:498
void newton(std::valarray< double > const &Dij, GradientProjection *gp, std::valarray< double > &coords, std::valarray< double > const &startCoords)
Definition cola.cpp:214
std::valarray< double > startX
Definition cola.h:472
vpsc::Rectangles boundingBoxes
Definition cola.h:463
void straighten(std::vector< straightener::Edge * > &, vpsc::Dim)
Definition cola.cpp:472
void setUnsatisfiableConstraintInfo(UnsatisfiableConstraintInfos *unsatisfiableX, UnsatisfiableConstraintInfos *unsatisfiableY)
Register to receive information about unsatisfiable constraints.
Definition cola.h:335
UnsatisfiableConstraintInfos * unsatisfiableX
Definition cola.h:499
NonOverlapConstraintsMode avoidOverlaps
Definition cola.h:500
void setStickyNodes(const double stickyWeight, std::valarray< double > const &startX, std::valarray< double > const &startY)
Sticky nodes causes nodes to spring back to (startX,startY) when unconstrained.
Definition cola.cpp:157
std::valarray< double > Y
Definition cola.h:469
void moveBoundingBoxes()
Update position of bounding boxes.
Definition cola.h:401
double compute_stress(std::valarray< double > const &Dij)
Definition cola.cpp:297
void setScaling(bool scaling)
Scaling speeds up the solver by conditioning the quadratic terms matrix.
Definition cola.h:352
GradientProjection * gpX
Definition cola.h:497
void setStraightenEdges(std::vector< straightener::Edge * > *straightenEdges, double bendWeight=0.01, double potBendWeight=0.1, bool xSkipping=true)
For the specified edges (with routings), generate dummy vars and constraints to try and straighten th...
Definition cola.h:385
std::vector< straightener::Edge * > * straightenEdges
Definition cola.h:501
void setConstraints(cola::CompoundConstraints *ccs)
Specify a set of compound constraints to apply to the layout.
Definition cola.h:306
UnsatisfiableConstraintInfos * unsatisfiableY
Definition cola.h:499
std::valarray< double > startY
Definition cola.h:473
void run(bool x=true, bool y=true)
Implements the main layout loop, taking descent steps until stress is no-longer significantly reduced...
Definition cola.cpp:319
void setNonOverlappingClusters()
Add constraints to prevent clusters overlapping.
Definition cola.h:375
void majorize(std::valarray< double > const &Dij, GradientProjection *gp, std::valarray< double > &coords, std::valarray< double > const &startCoords)
Definition cola.cpp:174
A Lock specifies a required position for a node.
Definition cola.h:78
Lock()
Definition cola.h:80
unsigned id
Definition cola.h:97
double x
Definition cola.h:98
Lock(unsigned id, double X, double Y)
Constructs a Lock object for a given node and position.
Definition cola.h:88
double y
Definition cola.h:99
unsigned getID() const
Definition cola.h:90
double pos(vpsc::Dim dim) const
Definition cola.h:93
A default functor that is called before each iteration in the main loop of the ConstrainedFDLayout::r...
Definition cola.h:168
virtual bool operator()()
Definition cola.h:196
static Resizes __resizesNotUsed
Definition cola.h:202
PreIteration(Locks &locks=__locksNotUsed, Resizes &resizes=__resizesNotUsed)
Constructs a PreIteration object that handles locking and resizing of nodes.
Definition cola.h:179
Resizes & resizes
Definition cola.h:198
Locks & locks
Definition cola.h:197
static Locks __locksNotUsed
Definition cola.h:201
PreIteration(Resizes &resizes)
Definition cola.h:185
virtual ~PreIteration()
Definition cola.h:194
A Resize specifies a new required bounding box for a node.
Definition cola.h:107
const vpsc::Rectangle * getTarget() const
Definition cola.h:127
unsigned id
Definition cola.h:131
Resize(unsigned id, double x, double y, double w, double h)
Constructs a Resize object for a given node and it's new bounding box.
Definition cola.h:122
vpsc::Rectangle target
Definition cola.h:132
unsigned getID() const
Definition cola.h:124
Holds the cluster hierarchy specification for a diagram.
Definition cluster.h:173
A default functor that is called after each iteration of the layout algorithm.
Definition cola.h:216
TestConvergence(const double tol=1e-4, const unsigned maxiterations=100)
Definition cola.h:219
const double tolerance
Definition cola.h:249
unsigned iterations
Definition cola.h:251
const unsigned maxiterations
Definition cola.h:250
virtual ~TestConvergence()
Definition cola.h:223
virtual bool operator()(const double new_stress, std::valarray< double > &X, std::valarray< double > &Y)
Definition cola.h:226
Interface for writing COLA addons to handle topology preserving layout.
Definition cola.h:532
virtual double applyForcesAndConstraints(ConstrainedFDLayout *layout, const vpsc::Dim dim, std::valarray< double > &g, vpsc::Variables &vs, vpsc::Constraints &cs, std::valarray< double > &coords, DesiredPositionsInDim &des, double oldStress)
Definition cola.h:595
virtual bool useTopologySolver(void) const
Definition cola.h:573
virtual void freeAssociatedObjects(void)
Definition cola.h:547
virtual double computeStress(void) const
Definition cola.h:569
virtual TopologyAddonInterface * clone(void) const
Definition cola.h:542
virtual void makeFeasible(bool generateNonOverlapConstraints, vpsc::Rectangles &boundingBoxes, cola::RootCluster *clusterHierarchy)
Definition cola.h:577
virtual void computePathLengths(unsigned short **G)
Definition cola.h:565
virtual void moveTo(const vpsc::Dim dim, vpsc::Variables &vs, vpsc::Constraints &cs, std::valarray< double > &coords, cola::RootCluster *clusterHierarchy)
Definition cola.h:585
virtual void handleResizes(const Resizes &resizeList, unsigned n, std::valarray< double > &X, std::valarray< double > &Y, cola::CompoundConstraints &ccs, vpsc::Rectangles &boundingBoxes, cola::RootCluster *clusterHierarchy)
Definition cola.h:551
virtual ~TopologyAddonInterface()
Definition cola.h:538
A rectangle represents a fixed-size shape in the diagram that may be moved to prevent overlaps and sa...
Definition rectangle.h:78
const double w
Definition conic-4.cpp:19
vector< Edge > es
vector< vpsc::Rectangle * > rs
double c[8][4]
libcola: Force-directed network layout subject to separation constraints library.
Definition box.cpp:25
ProjectionResult solve(vpsc::Variables &vs, vpsc::Constraints &cs, vpsc::Rectangles &rs, unsigned debugLevel=0)
Constructs a solver and attempts to solve the passed constraints on the passed vars.
std::vector< cola::Lock > Locks
A vector of Lock objects.
Definition cola.h:102
std::vector< CompoundConstraint * > CompoundConstraints
A vector of pointers to CompoundConstraint objects.
void setVariableDesiredPositions(vpsc::Variables &vs, vpsc::Constraints &cs, const DesiredPositionsInDim &des, std::valarray< double > &coords)
Definition colafd.cpp:1022
void removeClusterOverlapFast(RootCluster &clusterHierarchy, vpsc::Rectangles &rs, Locks &locks)
Definition cola.cpp:679
std::vector< cola::Resize > Resizes
A vector of Resize objects.
Definition cola.h:135
void setupVarsAndConstraints(unsigned n, const CompoundConstraints &ccs, const vpsc::Dim dim, vpsc::Rectangles &boundingBoxes, RootCluster *clusterHierarchy, vpsc::Variables &vs, vpsc::Constraints &cs, std::valarray< double > &coords)
Definition colafd.cpp:957
ConstrainedMajorizationLayout * simpleCMLFactory(vpsc::Rectangles &rs, std::vector< Edge > const &es, RootCluster *clusterHierarchy, const double idealLength, bool useNeighbourStress)
Definition cola.cpp:685
NonOverlapConstraintsMode
Definition commondefs.h:47
@ Both
Definition commondefs.h:47
@ Horizontal
Definition commondefs.h:47
std::vector< UnsatisfiableConstraintInfo * > UnsatisfiableConstraintInfos
A vector of pointers to UnsatisfiableConstraintInfo objects.
std::vector< unsigned > NodeIndexes
A vector of node Indexes.
Definition cola.h:62
std::vector< cola::DesiredPosition > DesiredPositions
Definition cola.h:147
std::vector< NodeIndexes > ListOfNodeIndexes
A vector of NodeIndexes.
Definition cola.h:65
void dijkstra(const unsigned s, const unsigned n, double *d, const std::vector< cola::Edge > &es, const std::valarray< double > &eLengths)
std::pair< unsigned, unsigned > Edge
Edges are simply a pair of indices to entries in the Node vector.
Definition cola.h:68
ProjectionResult projectOntoCCs(vpsc::Dim dim, vpsc::Rectangles &rs, cola::CompoundConstraints ccs, bool preventOverlaps, int accept=0, unsigned debugLevel=0)
Attempt to do a projection onto a vector of cola CompoundConstraints.
std::pair< unsigned, double > DesiredPositionInDim
Definition cola.h:154
std::vector< double > EdgeLengths
EdgeLengths is a vector of ideal lengths for edges corresponding to edges in the edge list.
Definition cola.h:72
std::vector< DesiredPositionInDim > DesiredPositionsInDim
Definition cola.h:155
Definition cola.h:44
libvpsc: Variable Placement with Separation Constraints quadratic program solver library.
Dim
Indicates the x- or y-dimension.
Definition rectangle.h:41
@ HORIZONTAL
The x-dimension (0).
Definition rectangle.h:43
std::vector< vpsc::Variable * > Variables
A vector of pointers to Variable objects.
std::vector< vpsc::Constraint * > Constraints
A vector of pointers to Constraint objects.
std::vector< Rectangle * > Rectangles
A vector of pointers to Rectangle objects.
Definition rectangle.h:246
std::string unsatinfo
Definition cola.h:927