Open
Graph Drawing
Framework

 v. 2023.09 (Elderberry)
 

FMMMLayout.h
Go to the documentation of this file.
1 
32 #pragma once
33 
34 #include <ogdf/basic/Graph.h>
36 #include <ogdf/basic/geometry.h>
40 
41 namespace ogdf {
42 
236 
237 public:
239  FMMMLayout();
240 
241  // destructor
242  virtual ~FMMMLayout() { }
243 
249  virtual void call(GraphAttributes& GA) override;
251 
255  void call(ClusterGraphAttributes& GA);
256 
258 
263  void call(GraphAttributes& GA, //graph and layout
264  const EdgeArray<double>& edgeLength); //factor for desired edge lengths
265 
267 
271  void call(GraphAttributes& GA, char* ps_file);
272 
274 
280  void call(GraphAttributes& GA, //graph and layout
281  const EdgeArray<double>& edgeLength, //factor for desired edge lengths
282  char* ps_file);
283 
289  double getCpuTime() { return time_total; }
291 
298 
302  void resetOptions();
303 
305 
319  bool useHighLevelOptions() const { return m_useHighLevelOptions; }
320 
322  void useHighLevelOptions(bool uho) { m_useHighLevelOptions = uho; }
323 
325  void setSingleLevel(bool b) { m_singleLevel = b; }
326 
328 
331  FMMMOptions::PageFormatType pageFormat() const { return m_pageFormat; }
332 
334  void pageFormat(FMMMOptions::PageFormatType t) { m_pageFormat = t; }
335 
337  double unitEdgeLength() const { return m_unitEdgeLength; }
338 
340  void unitEdgeLength(double x) { m_unitEdgeLength = ((x > 0.0) ? x : 1); }
341 
343 
351  bool newInitialPlacement() const { return m_newInitialPlacement; }
352 
354  void newInitialPlacement(bool nip) { m_newInitialPlacement = nip; }
355 
357 
361  FMMMOptions::QualityVsSpeed qualityVersusSpeed() const { return m_qualityVersusSpeed; }
362 
364  void qualityVersusSpeed(FMMMOptions::QualityVsSpeed qvs) { m_qualityVersusSpeed = qvs; }
365 
373  void randSeed(int p) { m_randSeed = ((0 <= p) ? p : 1); }
375 
377  int randSeed() const { return m_randSeed; }
378 
381  return m_edgeLengthMeasurement;
382  }
383 
386  m_edgeLengthMeasurement = elm;
387  }
388 
390  FMMMOptions::AllowedPositions allowedPositions() const { return m_allowedPositions; }
391 
393  void allowedPositions(FMMMOptions::AllowedPositions ap) { m_allowedPositions = ap; }
394 
396 
399  int maxIntPosExponent() const { return m_maxIntPosExponent; }
400 
402  void maxIntPosExponent(int e) { m_maxIntPosExponent = (((e >= 31) && (e <= 51)) ? e : 31); }
403 
409 
413  double pageRatio() const { return m_pageRatio; }
414 
416  void pageRatio(double r) { m_pageRatio = ((r > 0) ? r : 1); }
417 
419 
423  int stepsForRotatingComponents() const { return m_stepsForRotatingComponents; }
424 
426  void stepsForRotatingComponents(int n) { m_stepsForRotatingComponents = ((0 <= n) ? n : 0); }
427 
429  FMMMOptions::TipOver tipOverCCs() const { return m_tipOverCCs; }
430 
432  void tipOverCCs(FMMMOptions::TipOver to) { m_tipOverCCs = to; }
433 
435  double minDistCC() const { return m_minDistCC; }
436 
438  void minDistCC(double x) { m_minDistCC = ((x > 0) ? x : 1); }
439 
441  FMMMOptions::PreSort presortCCs() const { return m_presortCCs; }
442 
444  void presortCCs(FMMMOptions::PreSort ps) { m_presortCCs = ps; }
445 
451 
457  int minGraphSize() const { return m_minGraphSize; }
458 
460  void minGraphSize(int n) { m_minGraphSize = ((n >= 2) ? n : 2); }
461 
463  FMMMOptions::GalaxyChoice galaxyChoice() const { return m_galaxyChoice; }
464 
466  void galaxyChoice(FMMMOptions::GalaxyChoice gc) { m_galaxyChoice = gc; }
467 
469 
474  int randomTries() const { return m_randomTries; }
475 
477  void randomTries(int n) { m_randomTries = ((n >= 1) ? n : 1); }
478 
480  FMMMOptions::MaxIterChange maxIterChange() const { return m_maxIterChange; }
481 
483  void maxIterChange(FMMMOptions::MaxIterChange mic) { m_maxIterChange = mic; }
484 
486 
491  int maxIterFactor() const { return m_maxIterFactor; }
492 
494  void maxIterFactor(int f) { m_maxIterFactor = ((f >= 1) ? f : 1); }
495 
498  return m_initialPlacementMult;
499  }
500 
503  m_initialPlacementMult = ipm;
504  }
505 
511  FMMMOptions::ForceModel forceModel() const { return m_forceModel; }
513 
515  void forceModel(FMMMOptions::ForceModel fm) { m_forceModel = fm; }
516 
518  double springStrength() const { return m_springStrength; }
519 
521  void springStrength(double x) { m_springStrength = ((x > 0) ? x : 1); }
522 
524  double repForcesStrength() const { return m_repForcesStrength; }
525 
527  void repForcesStrength(double x) { m_repForcesStrength = ((x > 0) ? x : 1); }
528 
531  return m_repulsiveForcesCalculation;
532  }
533 
536  m_repulsiveForcesCalculation = rfc;
537  }
538 
540  FMMMOptions::StopCriterion stopCriterion() const { return m_stopCriterion; }
541 
543  void stopCriterion(FMMMOptions::StopCriterion rsc) { m_stopCriterion = rsc; }
544 
546 
550  double threshold() const { return m_threshold; }
551 
553  void threshold(double x) { m_threshold = ((x > 0) ? x : 0.1); }
554 
556  int fixedIterations() const { return m_fixedIterations; }
557 
559  void fixedIterations(int n) { m_fixedIterations = ((n >= 1) ? n : 1); }
560 
562  double forceScalingFactor() const { return m_forceScalingFactor; }
563 
565  void forceScalingFactor(double f) { m_forceScalingFactor = ((f > 0) ? f : 1); }
566 
568 
572  bool coolTemperature() const { return m_coolTemperature; }
573 
575  void coolTemperature(bool b) { m_coolTemperature = b; }
576 
578 
582  double coolValue() const { return m_coolValue; }
583 
585  void coolValue(double x) { m_coolValue = (((x > 0) && (x <= 1)) ? x : 0.99); }
586 
589  return m_initialPlacementForces;
590  }
591 
594  m_initialPlacementForces = ipf;
595  }
596 
602 
607  bool resizeDrawing() const { return m_resizeDrawing; }
608 
610  void resizeDrawing(bool b) { m_resizeDrawing = b; }
611 
613 
617  double resizingScalar() const { return m_resizingScalar; }
618 
620  void resizingScalar(double s) { m_resizingScalar = ((s > 0) ? s : 1); }
621 
623  int fineTuningIterations() const { return m_fineTuningIterations; }
624 
626  void fineTuningIterations(int n) { m_fineTuningIterations = ((n >= 0) ? n : 0); }
627 
629 
633  double fineTuneScalar() const { return m_fineTuneScalar; }
634 
636  void fineTuneScalar(double s) { m_fineTuneScalar = ((s >= 0) ? s : 1); }
637 
639 
644  bool adjustPostRepStrengthDynamically() const { return m_adjustPostRepStrengthDynamically; }
645 
647  void adjustPostRepStrengthDynamically(bool b) { m_adjustPostRepStrengthDynamically = b; }
648 
650  double postSpringStrength() const { return m_postSpringStrength; }
651 
653  void postSpringStrength(double x) { m_postSpringStrength = ((x > 0) ? x : 1); }
654 
656  double postStrengthOfRepForces() const { return m_postStrengthOfRepForces; }
657 
659  void postStrengthOfRepForces(double x) { m_postStrengthOfRepForces = ((x > 0) ? x : 1); }
660 
666 
671  int frGridQuotient() const { return m_frGridQuotient; }
672 
674  void frGridQuotient(int p) { m_frGridQuotient = ((0 <= p) ? p : 2); }
675 
677  FMMMOptions::ReducedTreeConstruction nmTreeConstruction() const { return m_NMTreeConstruction; }
678 
681  m_NMTreeConstruction = rtc;
682  }
683 
685  FMMMOptions::SmallestCellFinding nmSmallCell() const { return m_NMSmallCell; }
686 
688  void nmSmallCell(FMMMOptions::SmallestCellFinding scf) { m_NMSmallCell = scf; }
689 
691 
695  int nmParticlesInLeaves() const { return m_NMParticlesInLeaves; }
696 
698  void nmParticlesInLeaves(int n) { m_NMParticlesInLeaves = ((n >= 1) ? n : 1); }
699 
701  int nmPrecision() const { return m_NMPrecision; }
702 
704  void nmPrecision(int p) { m_NMPrecision = ((p >= 1) ? p : 1); }
705 
707 
708 private:
709  //high level options
715 
716  //low level options
717  //general options
722 
723  //options for divide et impera step
724  double m_pageRatio;
727  double m_minDistCC;
729 
730  //options for multilevel step
735 
741 
744 
745  //options for force calculation step
751  double m_threshold;
755  double m_coolValue;
757 
758  //options for postprocessing step
766 
767  //options for repulsive force approximation methods
774 
775  //other variables
777  double cool_factor;
779  double boxlength;
783  double time_total;
784 
787 
790 
792  void call_DIVIDE_ET_IMPERA_step(Graph& G, NodeArray<NodeAttributes>& A,
794 
796  void call_MULTILEVEL_step_for_subGraph(Graph& G, NodeArray<NodeAttributes>& A,
798 
800  bool running(int iter, int max_mult_iter, double actforcevectorlength);
801 
803 
808  void call_FORCE_CALCULATION_step(Graph& G, NodeArray<NodeAttributes>& A,
809  EdgeArray<EdgeAttributes>& E, int act_level, int max_level);
810 
812  void call_POSTPROCESSING_step(Graph& G, NodeArray<NodeAttributes>& A,
814  NodeArray<DPoint>& F_rep, NodeArray<DPoint>& last_node_movement);
815 
819 
821  void update_low_level_options_due_to_high_level_options_settings();
822 
824  void import_NodeAttributes(const Graph& G, GraphAttributes& GA, NodeArray<NodeAttributes>& A);
825 
827  void import_EdgeAttributes(const Graph& G, const EdgeArray<double>& edgeLength,
829 
831  void init_ind_ideal_edgelength(const Graph& G, NodeArray<NodeAttributes>& A,
833 
835  void set_radii(const Graph& G, NodeArray<NodeAttributes>& A);
836 
838  void export_NodeAttributes(Graph& G_reduced, NodeArray<NodeAttributes>& A_reduced,
839  GraphAttributes& GA);
840 
842 
847  void make_simple_loopfree(const Graph& G, NodeArray<NodeAttributes>& A,
849  EdgeArray<EdgeAttributes>& E_reduced);
850 
852 
856  void delete_parallel_edges(const Graph& G, EdgeArray<EdgeAttributes>& E, Graph& G_reduced,
857  List<edge>& S, EdgeArray<double>& new_edgelength);
858 
860 
863  void update_edgelength(List<edge>& S, EdgeArray<double>& new_edgelength,
864  EdgeArray<EdgeAttributes>& E_reduced);
865 
867 
870  double get_post_rep_force_strength(int n) { return min(0.2, 400.0 / double(n)); }
871 
877  void adjust_positions(const Graph& G, NodeArray<NodeAttributes>& A);
878 
880  void create_postscript_drawing(GraphAttributes& GA, char* ps_file);
881 
885 
887 
891  void create_maximum_connected_subGraphs(Graph& G, NodeArray<NodeAttributes>& A,
893  EdgeArray<EdgeAttributes> E_sub[], NodeArray<int>& component);
894 
896 
900  void pack_subGraph_drawings(NodeArray<NodeAttributes>& A, Graph G_sub[],
901  NodeArray<NodeAttributes> A_sub[]);
902 
904  void calculate_bounding_rectangles_of_components(List<Rectangle>& R, Graph G_sub[],
905  NodeArray<NodeAttributes> A_sub[]);
906 
908  Rectangle calculate_bounding_rectangle(Graph& G, NodeArray<NodeAttributes>& A,
909  int componenet_index);
910 
917  void rotate_components_and_calculate_bounding_rectangles(List<Rectangle>& R, Graph G_sub[],
918  NodeArray<NodeAttributes> A_sub[]);
919 
924  double calculate_area(double width, double height, int comp_nr) {
925  double scaling = 1.0;
926  if (comp_nr == 1) { //calculate aspect ratio area of the rectangle
927  OGDF_ASSERT(height != 0.0);
928  double ratio = width / height;
929  if (ratio < pageRatio()) { //scale width
930  OGDF_ASSERT(ratio != 0.0);
931  scaling = pageRatio() / ratio;
932  } else { //scale height
933  OGDF_ASSERT(pageRatio() != 0.0);
934  scaling = ratio / pageRatio();
935  }
936  }
937  return width * height * scaling;
938  }
939 
946  void export_node_positions(NodeArray<NodeAttributes>& A, List<Rectangle>& R, Graph G_sub[],
947  NodeArray<NodeAttributes> A_sub[]);
948 
951  EdgeArray<EdgeAttributes> E_sub[]) {
952  delete[] G_sub;
953  delete[] A_sub;
954  delete[] E_sub;
955  }
956 
960 
966  int get_max_mult_iter(int act_level, int max_level, int node_nr);
967 
971 
973  void calculate_forces(Graph& G, NodeArray<NodeAttributes>& A, EdgeArray<EdgeAttributes>& E,
975  NodeArray<DPoint>& last_node_movement, int iter, int fine_tuning_step);
976 
978  void init_boxlength_and_cornercoordinate(Graph& G, NodeArray<NodeAttributes>& A);
979 
981  void create_initial_placement(Graph& G, NodeArray<NodeAttributes>& A);
982 
984  void create_initial_placement_uniform_grid(const Graph& G, NodeArray<NodeAttributes>& A);
985 
987  void create_initial_placement_random(const Graph& G, NodeArray<NodeAttributes>& A);
988 
990  void init_F(Graph& G, NodeArray<DPoint>& F);
991 
992 
994  void make_initialisations_for_rep_calc_classes(
995  Graph& G/*,
996  NodeArray<NodeAttributes> &A,
997  NodeArray<DPoint>& F_rep*/);
998 
1001  switch (repulsiveForcesCalculation()) {
1003  FR.calculate_exact_repulsive_forces(G, A, F_rep);
1004  break;
1006  FR.calculate_approx_repulsive_forces(G, A, F_rep);
1007  break;
1009  NM.calculate_repulsive_forces(G, A, F_rep);
1010  break;
1011  }
1012  }
1013 
1016  if (repulsiveForcesCalculation() == FMMMOptions::RepulsiveForcesMethod::NMM) {
1017  NM.deallocate_memory();
1018  }
1019  }
1020 
1022  void calculate_attractive_forces(Graph& G, NodeArray<NodeAttributes>& A,
1024 
1026  double f_attr_scalar(double d, double ind_ideal_edge_length);
1027 
1029  void add_attr_rep_forces(Graph& G, NodeArray<DPoint>& F_attr, NodeArray<DPoint>& F_rep,
1030  NodeArray<DPoint>& F, int iter, int fine_tuning_step);
1031 
1034 
1036 
1039  void update_boxlength_and_cornercoordinate(Graph& G, NodeArray<NodeAttributes>& A);
1040 
1042  double max_radius(int iter) { return (iter == 1) ? boxlength / 1000 : boxlength / 5; }
1043 
1045  void set_average_ideal_edgelength(Graph& G, EdgeArray<EdgeAttributes>& E);
1046 
1051  double get_average_forcevector_length(Graph& G, NodeArray<DPoint>& F);
1052 
1057  void prevent_oscillations(Graph& G, NodeArray<DPoint>& F, NodeArray<DPoint>& last_node_movement,
1058  int iter);
1059 
1061  void init_last_node_movement(Graph& G, NodeArray<DPoint>& F,
1062  NodeArray<DPoint>& last_node_movement);
1063 
1068  void adapt_drawing_to_ideal_average_edgelength(Graph& G, NodeArray<NodeAttributes>& A,
1070 
1076  double x_min = down_left_corner.m_x;
1077  double x_max = down_left_corner.m_x + boxlength;
1078  double y_min = down_left_corner.m_y;
1079  double y_max = down_left_corner.m_y + boxlength;
1080  if (force.m_x < x_min) {
1081  force.m_x = x_min;
1082  } else if (force.m_x > x_max) {
1083  force.m_x = x_max;
1084  }
1085  if (force.m_y < y_min) {
1086  force.m_y = y_min;
1087  } else if (force.m_y > y_max) {
1088  force.m_y = y_max;
1089  }
1090  }
1091 
1095 
1097  void init_time() { time_total = 0; }
1098 
1100 };
1101 
1102 }
ogdf::FMMMLayout::stopCriterion
FMMMOptions::StopCriterion stopCriterion() const
Returns the stop criterion.
Definition: FMMMLayout.h:540
ogdf::FMMMLayout::fineTuneScalar
void fineTuneScalar(double s)
Sets the option fineTuneScalar to s.
Definition: FMMMLayout.h:636
ogdf::FMMMLayout::down_left_corner
DPoint down_left_corner
Holds down left corner of the comput.
Definition: FMMMLayout.h:781
ogdf::FMMMLayout::nmSmallCell
void nmSmallCell(FMMMOptions::SmallestCellFinding scf)
Sets the option nmSmallCell to scf.
Definition: FMMMLayout.h:688
ogdf::FMMMLayout::m_postStrengthOfRepForces
double m_postStrengthOfRepForces
The strength of repulsive forces during postprocessing.
Definition: FMMMLayout.h:765
ogdf::FMMMLayout::m_tipOverCCs
FMMMOptions::TipOver m_tipOverCCs
Option for tip-over of connected components.
Definition: FMMMLayout.h:726
ogdf
The namespace for all OGDF objects.
Definition: AugmentationModule.h:36
ogdf::GraphAttributes
Stores additional attributes of a graph (like layout information).
Definition: GraphAttributes.h:66
ogdf::FMMMLayout::m_NMParticlesInLeaves
int m_NMParticlesInLeaves
The maximal number of particles in a leaf.
Definition: FMMMLayout.h:772
ogdf::FMMMLayout::maxIterFactor
int maxIterFactor() const
Returns the current setting of option maxIterFactor.
Definition: FMMMLayout.h:491
ogdf::FMMMLayout::resizeDrawing
void resizeDrawing(bool b)
Sets the option resizeDrawing to b.
Definition: FMMMLayout.h:610
ogdf::FMMMLayout::m_threshold
double m_threshold
The threshold for the stop criterion.
Definition: FMMMLayout.h:751
ogdf::FMMMLayout::m_unitEdgeLength
double m_unitEdgeLength
The unit edge length.
Definition: FMMMLayout.h:712
ogdf::FMMMLayout::m_maxIntPosExponent
int m_maxIntPosExponent
The option for the used exponent.
Definition: FMMMLayout.h:721
ogdf::FMMMLayout::fixedIterations
int fixedIterations() const
Returns the fixed number of iterations for the stop criterion.
Definition: FMMMLayout.h:556
ogdf::FMMMLayout::edgeLengthMeasurement
FMMMOptions::EdgeLengthMeasurement edgeLengthMeasurement() const
Returns the current setting of option edgeLengthMeasurement.
Definition: FMMMLayout.h:380
Graph.h
Includes declaration of graph class.
ogdf::GenericPoint< double >
ogdf::FMMMOptions::RepulsiveForcesMethod::NMM
@ NMM
Calculation as for new multipole method (fast and accurate).
ogdf::FMMMLayout::galaxyChoice
void galaxyChoice(FMMMOptions::GalaxyChoice gc)
Sets the option galaxyChoice to gc.
Definition: FMMMLayout.h:466
ogdf::FMMMLayout::setSingleLevel
void setSingleLevel(bool b)
Sets single level option, no multilevel hierarchy is created if b == true.
Definition: FMMMLayout.h:325
ogdf::FMMMLayout::repulsiveForcesCalculation
FMMMOptions::RepulsiveForcesMethod repulsiveForcesCalculation() const
Returns the current setting of option repulsiveForcesCalculation.
Definition: FMMMLayout.h:530
ogdf::FMMMLayout::m_minDistCC
double m_minDistCC
The separation between connected components.
Definition: FMMMLayout.h:727
geometry.h
Declaration of classes GenericPoint, GenericPolyline, GenericLine, GenericSegment,...
ogdf::energybased::fmmm::NewMultipoleMethod::deallocate_memory
void deallocate_memory()
Dynamically allocated memory is freed here.
ogdf::FMMMLayout::qualityVersusSpeed
FMMMOptions::QualityVsSpeed qualityVersusSpeed() const
Returns the current setting of option qualityVersusSpeed.
Definition: FMMMLayout.h:361
ogdf::FMMMLayout::stepsForRotatingComponents
int stepsForRotatingComponents() const
Returns the current setting of option stepsForRotatingComponents.
Definition: FMMMLayout.h:423
OGDF_ASSERT
#define OGDF_ASSERT(expr)
Assert condition expr. See doc/build.md for more information.
Definition: basic.h:54
ogdf::FMMMLayout::nmPrecision
void nmPrecision(int p)
Sets the precision for the multipole expansions to p.
Definition: FMMMLayout.h:704
ogdf::FMMMOptions::PreSort
PreSort
Specifies how connected components are sorted before the packing algorithm is applied.
Definition: FMMMOptions.h:80
ogdf::FMMMLayout::m_edgeLengthMeasurement
FMMMOptions::EdgeLengthMeasurement m_edgeLengthMeasurement
The option for edge length measurement.
Definition: FMMMLayout.h:719
ogdf::FMMMLayout::m_stopCriterion
FMMMOptions::StopCriterion m_stopCriterion
The stop criterion.
Definition: FMMMLayout.h:750
ogdf::FMMMLayout
The fast multipole multilevel layout algorithm.
Definition: FMMMLayout.h:232
ogdf::FMMMLayout::max_integer_position
double max_integer_position
The maximum value for an integer position.
Definition: FMMMLayout.h:776
ogdf::FMMMLayout::adjustPostRepStrengthDynamically
bool adjustPostRepStrengthDynamically() const
Returns the current setting of option adjustPostRepStrengthDynamically.
Definition: FMMMLayout.h:644
ogdf::FMMMLayout::m_frGridQuotient
int m_frGridQuotient
The grid quotient.
Definition: FMMMLayout.h:768
ogdf::FMMMLayout::allowedPositions
FMMMOptions::AllowedPositions allowedPositions() const
Returns the current setting of option allowedPositions.
Definition: FMMMLayout.h:390
ogdf::FMMMLayout::postStrengthOfRepForces
double postStrengthOfRepForces() const
Returns the strength of the repulsive forces in the postprocessing step.
Definition: FMMMLayout.h:656
ogdf::GenericPoint::m_y
T m_y
The y-coordinate.
Definition: geometry.h:80
ogdf::FMMMLayout::initialPlacementForces
void initialPlacementForces(FMMMOptions::InitialPlacementForces ipf)
Sets the option initialPlacementForces to ipf.
Definition: FMMMLayout.h:593
ogdf::FMMMLayout::frGridQuotient
int frGridQuotient() const
Returns the current setting of option frGridQuotient.
Definition: FMMMLayout.h:671
ogdf::FMMMOptions::GalaxyChoice
GalaxyChoice
Specifies how sun nodes of galaxies are selected.
Definition: FMMMOptions.h:87
ogdf::FMMMLayout::m_allowedPositions
FMMMOptions::AllowedPositions m_allowedPositions
The option for allowed positions.
Definition: FMMMLayout.h:720
ogdf::FMMMLayout::forceScalingFactor
void forceScalingFactor(double f)
Sets the scaling factor for the forces to f.
Definition: FMMMLayout.h:565
ogdf::whaType::A
@ A
ogdf::FMMMLayout::maxIntPosExponent
int maxIntPosExponent() const
Returns the current setting of option maxIntPosExponent.
Definition: FMMMLayout.h:399
ogdf::FMMMLayout::resizingScalar
double resizingScalar() const
Returns the current setting of option resizingScalar.
Definition: FMMMLayout.h:617
ogdf::FMMMLayout::fineTuningIterations
int fineTuningIterations() const
Returns the number of iterations for fine tuning.
Definition: FMMMLayout.h:623
ogdf::FMMMLayout::springStrength
void springStrength(double x)
Sets the strength of the springs to x.
Definition: FMMMLayout.h:521
ogdf::FMMMLayout::threshold
double threshold() const
Returns the threshold for the stop criterion.
Definition: FMMMLayout.h:550
ogdf::FMMMLayout::average_ideal_edgelength
double average_ideal_edgelength
Measured from center to center.
Definition: FMMMLayout.h:778
ogdf::FMMMLayout::randomTries
void randomTries(int n)
Sets the option randomTries to n.
Definition: FMMMLayout.h:477
ogdf::FMMMLayout::m_maxIterChange
FMMMOptions::MaxIterChange m_maxIterChange
The option for how to change MaxIterations. If maxIterChange != micConstant, the iterations are decre...
Definition: FMMMLayout.h:740
ogdf::energybased::fmmm::NewMultipoleMethod
Definition: NewMultipoleMethod.h:48
ogdf::FMMMLayout::coolTemperature
void coolTemperature(bool b)
Sets the option coolTemperature to b.
Definition: FMMMLayout.h:575
ogdf::energybased::fmmm::NodeAttributes
helping data structure that stores the graphical attributes of a node that are needed for the force-d...
Definition: NodeAttributes.h:44
LayoutModule.h
Declaration of interface for layout algorithms (class LayoutModule)
ogdf::FMMMLayout::m_presortCCs
FMMMOptions::PreSort m_presortCCs
The option for presorting connected components.
Definition: FMMMLayout.h:728
ogdf::FMMMLayout::springStrength
double springStrength() const
Returns the strength of the springs.
Definition: FMMMLayout.h:518
ogdf::FMMMLayout::resizeDrawing
bool resizeDrawing() const
Returns the current setting of option resizeDrawing.
Definition: FMMMLayout.h:607
ogdf::FMMMLayout::nmParticlesInLeaves
int nmParticlesInLeaves() const
Returns the current setting of option nmParticlesInLeaves.
Definition: FMMMLayout.h:695
ogdf::FMMMLayout::nmPrecision
int nmPrecision() const
Returns the precision p for the p-term multipole expansions.
Definition: FMMMLayout.h:701
ogdf::FMMMLayout::boxlength
double boxlength
Holds the length of the quadratic comput.
Definition: FMMMLayout.h:779
ogdf::FMMMLayout::frGridQuotient
void frGridQuotient(int p)
Sets the option frGridQuotient to p.
Definition: FMMMLayout.h:674
ogdf::FMMMLayout::m_initialPlacementMult
FMMMOptions::InitialPlacementMult m_initialPlacementMult
The option for creating initial placement.
Definition: FMMMLayout.h:743
ogdf::FMMMOptions::InitialPlacementForces
InitialPlacementForces
Specifies how the initial placement is done.
Definition: FMMMOptions.h:129
ogdf::FMMMLayout::minGraphSize
int minGraphSize() const
Returns the current setting of option minGraphSize.
Definition: FMMMLayout.h:457
ogdf::FMMMLayout::pageRatio
void pageRatio(double r)
Sets the option pageRatio to r.
Definition: FMMMLayout.h:416
ogdf::FMMMLayout::m_fineTuningIterations
int m_fineTuningIterations
The number of iterations for fine tuning.
Definition: FMMMLayout.h:761
ogdf::FMMMLayout::m_useHighLevelOptions
bool m_useHighLevelOptions
The option for using high-level options.
Definition: FMMMLayout.h:710
NewMultipoleMethod.h
Declaration of class NewMultipoleMethod (NMM).
ogdf::FMMMOptions::InitialPlacementMult
InitialPlacementMult
Specifies how the initial placement is generated.
Definition: FMMMOptions.h:102
ogdf::FMMMLayout::pageRatio
double pageRatio() const
Returns the current setting of option pageRatio.
Definition: FMMMLayout.h:413
ogdf::FMMMLayout::initialPlacementMult
FMMMOptions::InitialPlacementMult initialPlacementMult() const
Returns the current setting of option initialPlacementMult.
Definition: FMMMLayout.h:497
ogdf::FMMMLayout::m_coolValue
double m_coolValue
The value by which forces are decreased.
Definition: FMMMLayout.h:755
ogdf::FMMMLayout::fixedIterations
void fixedIterations(int n)
Sets the fixed number of iterations for the stop criterion to n.
Definition: FMMMLayout.h:559
ogdf::FMMMOptions::MaxIterChange
MaxIterChange
Specifies how MaxIterations is changed in subsequent multilevels.
Definition: FMMMOptions.h:95
ogdf::FMMMLayout::m_newInitialPlacement
bool m_newInitialPlacement
The option for new initial placement.
Definition: FMMMLayout.h:713
ogdf::FMMMLayout::unitEdgeLength
double unitEdgeLength() const
Returns the current setting of option unitEdgeLength.
Definition: FMMMLayout.h:337
ogdf::GenericPoint::m_x
T m_x
The x-coordinate.
Definition: geometry.h:79
ogdf::FMMMOptions::AllowedPositions
AllowedPositions
Specifies which positions for a node are allowed.
Definition: FMMMOptions.h:59
ogdf::FMMMLayout::pageFormat
void pageFormat(FMMMOptions::PageFormatType t)
Sets the option pageRatio to t.
Definition: FMMMLayout.h:334
ogdf::FMMMLayout::useHighLevelOptions
bool useHighLevelOptions() const
Returns the current setting of option useHighLevelOptions.
Definition: FMMMLayout.h:319
ogdf::FMMMLayout::minGraphSize
void minGraphSize(int n)
Sets the option minGraphSize to n.
Definition: FMMMLayout.h:460
ogdf::FMMMLayout::m_fixedIterations
int m_fixedIterations
The fixed number of iterations for the stop criterion.
Definition: FMMMLayout.h:752
ogdf::FMMMLayout::m_qualityVersusSpeed
FMMMOptions::QualityVsSpeed m_qualityVersusSpeed
The option for quality-vs-speed trade-off.
Definition: FMMMLayout.h:714
ogdf::ClusterGraphAttributes
Stores additional attributes of a clustered graph (like layout information).
Definition: ClusterGraphAttributes.h:46
ogdf::FMMMLayout::m_fineTuneScalar
double m_fineTuneScalar
Parameter for scaling forces during fine tuning.
Definition: FMMMLayout.h:762
ogdf::FMMMLayout::maxIterChange
FMMMOptions::MaxIterChange maxIterChange() const
Returns the current setting of option maxIterChange.
Definition: FMMMLayout.h:480
r
int r[]
Definition: hierarchical-ranking.cpp:8
ogdf::energybased::fmmm::FruchtermanReingold::calculate_approx_repulsive_forces
void calculate_approx_repulsive_forces(const Graph &G, NodeArray< NodeAttributes > &A, NodeArray< DPoint > &F_rep)
Grid approximation of rep.forces for each node.
ogdf::FMMMOptions::StopCriterion
StopCriterion
Specifies the stop criterion.
Definition: FMMMOptions.h:122
ogdf::energybased::fmmm::Rectangle
Helping data structure for packing rectangles; The width, height and the position of the down left co...
Definition: Rectangle.h:45
ogdf::FMMMLayout::m_randSeed
int m_randSeed
The random seed.
Definition: FMMMLayout.h:718
ogdf::FMMMLayout::calculate_repulsive_forces
void calculate_repulsive_forces(Graph &G, NodeArray< NodeAttributes > &A, NodeArray< DPoint > &F_rep)
Calculates repulsive forces for each node.
Definition: FMMMLayout.h:1000
ogdf::FMMMLayout::m_pageFormat
FMMMOptions::PageFormatType m_pageFormat
The option for the page format.
Definition: FMMMLayout.h:711
ogdf::FMMMLayout::m_randomTries
int m_randomTries
The number of random tries.
Definition: FMMMLayout.h:734
ogdf::FMMMLayout::m_initialPlacementForces
FMMMOptions::InitialPlacementForces m_initialPlacementForces
The option for how the initial placement is done.
Definition: FMMMLayout.h:756
ogdf::FMMMLayout::m_postSpringStrength
double m_postSpringStrength
The strength of springs during postprocessing.
Definition: FMMMLayout.h:764
ogdf::FMMMLayout::presortCCs
FMMMOptions::PreSort presortCCs() const
Returns the current setting of option presortCCs.
Definition: FMMMLayout.h:441
ogdf::FMMMLayout::initialPlacementForces
FMMMOptions::InitialPlacementForces initialPlacementForces() const
Returns the current setting of option initialPlacementForces.
Definition: FMMMLayout.h:588
ogdf::FMMMOptions::TipOver
TipOver
Specifies in which case it is allowed to tip over drawings of connected components.
Definition: FMMMOptions.h:73
ogdf::FMMMOptions::ForceModel
ForceModel
Specifies the force model.
Definition: FMMMOptions.h:108
ogdf::FMMMLayout::maxIterFactor
void maxIterFactor(int f)
Sets the option maxIterFactor to f.
Definition: FMMMLayout.h:494
ogdf::FMMMLayout::repForcesStrength
double repForcesStrength() const
Returns the strength of the repulsive forces.
Definition: FMMMLayout.h:524
ogdf::energybased::fmmm::FruchtermanReingold
Definition: FruchtermanReingold.h:44
ogdf::FMMMLayout::m_pageRatio
double m_pageRatio
The desired page ratio.
Definition: FMMMLayout.h:724
ogdf::FMMMLayout::maxIterChange
void maxIterChange(FMMMOptions::MaxIterChange mic)
Sets the option maxIterChange to mic.
Definition: FMMMLayout.h:483
ogdf::FMMMLayout::pageFormat
FMMMOptions::PageFormatType pageFormat() const
Returns the current setting of option pageFormat.
Definition: FMMMLayout.h:331
ogdf::FMMMLayout::restrict_force_to_comp_box
void restrict_force_to_comp_box(DPoint &force)
The force is restricted to have values within the comp.
Definition: FMMMLayout.h:1075
ogdf::FMMMOptions::PageFormatType
PageFormatType
Possible page formats.
Definition: FMMMOptions.h:39
ogdf::FMMMLayout::radius
NodeArray< double > radius
Holds the radius of the surrounding circle for each node.
Definition: FMMMLayout.h:782
ogdf::FMMMLayout::allowedPositions
void allowedPositions(FMMMOptions::AllowedPositions ap)
Sets the option allowedPositions to ap.
Definition: FMMMLayout.h:393
ogdf::FMMMLayout::calculate_area
double calculate_area(double width, double height, int comp_nr)
Returns the area (aspect ratio area) of a rectangle with width w and height h if comp_nr > 1 ( comp_n...
Definition: FMMMLayout.h:924
ogdf::FMMMLayout::repulsiveForcesCalculation
void repulsiveForcesCalculation(FMMMOptions::RepulsiveForcesMethod rfc)
Sets the option repulsiveForcesCalculation to rfc.
Definition: FMMMLayout.h:535
ogdf::FMMMLayout::newInitialPlacement
void newInitialPlacement(bool nip)
Sets the option newInitialPlacement to nip.
Definition: FMMMLayout.h:354
ogdf::FMMMLayout::nmTreeConstruction
FMMMOptions::ReducedTreeConstruction nmTreeConstruction() const
Returns the current setting of option nmTreeConstruction.
Definition: FMMMLayout.h:677
ogdf::FMMMLayout::initialPlacementMult
void initialPlacementMult(FMMMOptions::InitialPlacementMult ipm)
Sets the option initialPlacementMult to ipm.
Definition: FMMMLayout.h:502
ogdf::FMMMLayout::nmTreeConstruction
void nmTreeConstruction(FMMMOptions::ReducedTreeConstruction rtc)
Sets the option nmTreeConstruction to rtc.
Definition: FMMMLayout.h:680
ogdf::FMMMLayout::postSpringStrength
void postSpringStrength(double x)
Sets the strength of the springs in the postprocessing step to x.
Definition: FMMMLayout.h:653
ogdf::FMMMLayout::presortCCs
void presortCCs(FMMMOptions::PreSort ps)
Sets the option presortCCs to ps.
Definition: FMMMLayout.h:444
ogdf::FMMMLayout::m_repForcesStrength
double m_repForcesStrength
The strength of repulsive forces.
Definition: FMMMLayout.h:748
ogdf::List< edge >
ogdf::internal::GraphRegisteredArray
RegisteredArray for nodes, edges and adjEntries of a graph.
Definition: Graph_d.h:651
ogdf::FMMMLayout::minDistCC
void minDistCC(double x)
Sets the minimal distance between connected components to x.
Definition: FMMMLayout.h:438
ogdf::energybased::fmmm::FruchtermanReingold::calculate_exact_repulsive_forces
void calculate_exact_repulsive_forces(const Graph &G, NodeArray< NodeAttributes > &A, NodeArray< DPoint > &F_rep)
Calculate exact rep. forces for each node.
ogdf::Graph
Data type for general directed graphs (adjacency list representation).
Definition: Graph_d.h:862
ogdf::FMMMOptions::EdgeLengthMeasurement
EdgeLengthMeasurement
Specifies how the length of an edge is measured.
Definition: FMMMOptions.h:53
ogdf::FMMMLayout::deallocate_memory_for_rep_calc_classes
void deallocate_memory_for_rep_calc_classes()
Deallocates dynamically allocated memory of the choosen rep. calculation class.
Definition: FMMMLayout.h:1015
ogdf::FMMMLayout::coolValue
void coolValue(double x)
Sets the option coolValue to x.
Definition: FMMMLayout.h:585
Rectangle.h
Declaration of class Rectangle.
ogdf::FMMMLayout::m_resizeDrawing
bool m_resizeDrawing
The option for resizing the drawing.
Definition: FMMMLayout.h:759
ogdf::FMMMLayout::m_repulsiveForcesCalculation
FMMMOptions::RepulsiveForcesMethod m_repulsiveForcesCalculation
Option for how to calculate repulsive forces.
Definition: FMMMLayout.h:749
ogdf::FMMMLayout::m_forceScalingFactor
double m_forceScalingFactor
The scaling factor for the forces.
Definition: FMMMLayout.h:753
ogdf::FMMMLayout::m_galaxyChoice
FMMMOptions::GalaxyChoice m_galaxyChoice
The selection of galaxy nodes.
Definition: FMMMLayout.h:733
ogdf::FMMMLayout::nmSmallCell
FMMMOptions::SmallestCellFinding nmSmallCell() const
Returns the current setting of option nmSmallCell.
Definition: FMMMLayout.h:685
ogdf::FMMMLayout::m_forceModel
FMMMOptions::ForceModel m_forceModel
The used force model.
Definition: FMMMLayout.h:746
ogdf::FMMMLayout::nmParticlesInLeaves
void nmParticlesInLeaves(int n)
Sets the option nmParticlesInLeaves to n.
Definition: FMMMLayout.h:698
ogdf::fast_multipole_embedder::move_nodes
double move_nodes(float *x, float *y, const uint32_t begin, const uint32_t end, const float *fx, const float *fy, const float t)
Definition: FMEKernel.h:77
ogdf::FMMMLayout::coolValue
double coolValue() const
Returns the current setting of option coolValue.
Definition: FMMMLayout.h:582
ogdf::FMMMLayout::delete_all_subGraphs
void delete_all_subGraphs(Graph G_sub[], NodeArray< NodeAttributes > A_sub[], EdgeArray< EdgeAttributes > E_sub[])
Frees dynamically allocated memory for the connected component subgraphs.
Definition: FMMMLayout.h:950
ogdf::energybased::fmmm::EdgeAttributes
helping data structure that stores the graphical attributes of an edge that are needed for the force-...
Definition: EdgeAttributes.h:43
ogdf::FMMMLayout::init_time
void init_time()
Sets time_total to zero.
Definition: FMMMLayout.h:1097
ogdf::FMMMLayout::threshold
void threshold(double x)
Sets the threshold for the stop criterion to x.
Definition: FMMMLayout.h:553
ogdf::FMMMLayout::stopCriterion
void stopCriterion(FMMMOptions::StopCriterion rsc)
Sets the stop criterion to rsc.
Definition: FMMMLayout.h:543
ogdf::FMMMLayout::number_of_components
int number_of_components
The number of components of the graph.
Definition: FMMMLayout.h:780
ogdf::FMMMLayout::postStrengthOfRepForces
void postStrengthOfRepForces(double x)
Sets the strength of the repulsive forces in the postprocessing step to x.
Definition: FMMMLayout.h:659
ogdf::FMMMLayout::m_NMTreeConstruction
FMMMOptions::ReducedTreeConstruction m_NMTreeConstruction
The option for how to construct reduced bucket quadtree.
Definition: FMMMLayout.h:770
ogdf::FMMMLayout::m_stepsForRotatingComponents
int m_stepsForRotatingComponents
The number of rotations.
Definition: FMMMLayout.h:725
ogdf::FMMMOptions::RepulsiveForcesMethod::Exact
@ Exact
Exact calculation (slow).
ogdf::FMMMLayout::newInitialPlacement
bool newInitialPlacement() const
Returns the current setting of option newInitialPlacement.
Definition: FMMMLayout.h:351
ogdf::FMMMLayout::m_adjustPostRepStrengthDynamically
bool m_adjustPostRepStrengthDynamically
The option adjustPostRepStrengthDynamically.
Definition: FMMMLayout.h:763
ogdf::FMMMLayout::FR
energybased::fmmm::FruchtermanReingold FR
Class for repulsive force calculation (Fruchterman, Reingold).
Definition: FMMMLayout.h:785
ogdf::FMMMLayout::get_post_rep_force_strength
double get_post_rep_force_strength(int n)
Returns the value for the strength of the repulsive forces.
Definition: FMMMLayout.h:870
ogdf::FMMMLayout::~FMMMLayout
virtual ~FMMMLayout()
Definition: FMMMLayout.h:242
ogdf::FMMMOptions::SmallestCellFinding
SmallestCellFinding
Specifies how to calculate the smallest quadratic cell that surrounds the particles of a node in the ...
Definition: FMMMOptions.h:144
ogdf::FMMMLayout::randSeed
int randSeed() const
Returns the seed of the random number generator.
Definition: FMMMLayout.h:377
ogdf::FMMMLayout::max_radius
double max_radius(int iter)
Describes the max. radius of a move in one time step, depending on the number of iterations.
Definition: FMMMLayout.h:1042
OGDF_EXPORT
#define OGDF_EXPORT
Specifies that a function or class is exported by the OGDF DLL.
Definition: config.h:101
ogdf::FMMMLayout::m_NMSmallCell
FMMMOptions::SmallestCellFinding m_NMSmallCell
The option for how to calculate smallest quadtratic cells.
Definition: FMMMLayout.h:771
ogdf::FMMMLayout::maxIntPosExponent
void maxIntPosExponent(int e)
Sets the option maxIntPosExponent to e.
Definition: FMMMLayout.h:402
ogdf::energybased::fmmm::NewMultipoleMethod::calculate_repulsive_forces
void calculate_repulsive_forces(const Graph &G, NodeArray< NodeAttributes > &A, NodeArray< DPoint > &F_rep)
Calculate rep. forces for each node.
ogdf::FMMMLayout::qualityVersusSpeed
void qualityVersusSpeed(FMMMOptions::QualityVsSpeed qvs)
Sets the option qualityVersusSpeed to qvs.
Definition: FMMMLayout.h:364
ogdf::FMMMLayout::repForcesStrength
void repForcesStrength(double x)
Sets the strength of the repulsive forces to x.
Definition: FMMMLayout.h:527
ogdf::FMMMLayout::time_total
double time_total
The runtime (=CPU-time) of the algorithm in seconds.
Definition: FMMMLayout.h:783
ogdf::FMMMLayout::m_springStrength
double m_springStrength
The strengths of springs.
Definition: FMMMLayout.h:747
ogdf::FMMMLayout::m_coolTemperature
bool m_coolTemperature
The option for how to scale forces.
Definition: FMMMLayout.h:754
ogdf::FMMMLayout::tipOverCCs
FMMMOptions::TipOver tipOverCCs() const
Returns the current setting of option tipOverCCs.
Definition: FMMMLayout.h:429
ogdf::FMMMOptions::RepulsiveForcesMethod
RepulsiveForcesMethod
Specifies how to calculate repulsive forces.
Definition: FMMMOptions.h:115
ogdf::FMMMLayout::cool_factor
double cool_factor
Needed for scaling the forces if coolTemperature is true.
Definition: FMMMLayout.h:777
ClusterGraphAttributes.h
Declares ClusterGraphAttributes, an extension of class GraphAttributes, to store clustergraph layout ...
ogdf::FMMMLayout::fineTuneScalar
double fineTuneScalar() const
Returns the curent setting of option fineTuneScalar.
Definition: FMMMLayout.h:633
ogdf::FMMMLayout::tipOverCCs
void tipOverCCs(FMMMOptions::TipOver to)
Sets the option tipOverCCs to to.
Definition: FMMMLayout.h:432
ogdf::FMMMLayout::postSpringStrength
double postSpringStrength() const
Returns the strength of the springs in the postprocessing step.
Definition: FMMMLayout.h:650
ogdf::FMMMLayout::unitEdgeLength
void unitEdgeLength(double x)
Sets the option unitEdgeLength to x.
Definition: FMMMLayout.h:340
ogdf::FMMMOptions::ReducedTreeConstruction
ReducedTreeConstruction
Specifies how the reduced bucket quadtree is constructed.
Definition: FMMMOptions.h:137
ogdf::FMMMLayout::randomTries
int randomTries() const
Returns the current setting of option randomTries.
Definition: FMMMLayout.h:474
ogdf::FMMMLayout::useHighLevelOptions
void useHighLevelOptions(bool uho)
Sets the option useHighLevelOptions to uho.
Definition: FMMMLayout.h:322
ogdf::FMMMLayout::coolTemperature
bool coolTemperature() const
Returns the current setting of option coolTemperature.
Definition: FMMMLayout.h:572
ogdf::FMMMOptions::QualityVsSpeed
QualityVsSpeed
Trade-off between run-time and quality.
Definition: FMMMOptions.h:46
ogdf::FMMMLayout::m_minGraphSize
int m_minGraphSize
The option for minimal graph size.
Definition: FMMMLayout.h:732
ogdf::FMMMLayout::edgeLengthMeasurement
void edgeLengthMeasurement(FMMMOptions::EdgeLengthMeasurement elm)
Sets the option edgeLengthMeasurement to elm.
Definition: FMMMLayout.h:385
ogdf::FMMMLayout::forceModel
void forceModel(FMMMOptions::ForceModel fm)
Sets the used force model to fm.
Definition: FMMMLayout.h:515
ogdf::FMMMLayout::minDistCC
double minDistCC() const
Returns the minimal distance between connected components.
Definition: FMMMLayout.h:435
ogdf::FMMMLayout::m_maxIterFactor
int m_maxIterFactor
The factor used for decreasing MaxIterations.
Definition: FMMMLayout.h:742
ogdf::FMMMLayout::m_resizingScalar
double m_resizingScalar
Parameter for resizing the drawing.
Definition: FMMMLayout.h:760
ogdf::FMMMLayout::forceScalingFactor
double forceScalingFactor() const
Returns the scaling factor for the forces.
Definition: FMMMLayout.h:562
ogdf::FMMMLayout::fineTuningIterations
void fineTuningIterations(int n)
Sets the number of iterations for fine tuning to n.
Definition: FMMMLayout.h:626
ogdf::FMMMLayout::m_singleLevel
bool m_singleLevel
Option for pure single level.
Definition: FMMMLayout.h:731
ogdf::FMMMLayout::m_NMPrecision
int m_NMPrecision
The precision for multipole expansions.
Definition: FMMMLayout.h:773
ogdf::FMMMLayout::galaxyChoice
FMMMOptions::GalaxyChoice galaxyChoice() const
Returns the current setting of option galaxyChoice.
Definition: FMMMLayout.h:463
ogdf::FMMMLayout::resizingScalar
void resizingScalar(double s)
Sets the option resizingScalar to s.
Definition: FMMMLayout.h:620
ogdf::FMMMLayout::stepsForRotatingComponents
void stepsForRotatingComponents(int n)
Sets the option stepsForRotatingComponents to n.
Definition: FMMMLayout.h:426
ogdf::FMMMLayout::adjustPostRepStrengthDynamically
void adjustPostRepStrengthDynamically(bool b)
Sets the option adjustPostRepStrengthDynamically to b.
Definition: FMMMLayout.h:647
ogdf::FMMMOptions::RepulsiveForcesMethod::GridApproximation
@ GridApproximation
Grid approximation (inaccurate).
ogdf::FMMMLayout::NM
energybased::fmmm::NewMultipoleMethod NM
Class for repulsive force calculation.
Definition: FMMMLayout.h:786
ogdf::internal::EdgeArrayBase2
RegisteredArray for edges of a graph, specialized for EdgeArray<edge>.
Definition: Graph_d.h:709
ogdf::LayoutModule
Interface of general layout algorithms.
Definition: LayoutModule.h:44