WaveBEM: Unsteady Nonlinear Potential Flow Solver for Ship-Wave Interaction.
occ_line_smoothing.h
Go to the documentation of this file.
1 #ifndef occ_line_smoothing_h
2 #define occ_line_smoothing_h
3 
4 #include <string>
5 #include <TopoDS_Shape.hxx>
6 #include <Geom_Plane.hxx>
7 #include <GeomLib_Tool.hxx>
8 #include <GeomAdaptor_Curve.hxx>
9 #include <gp_Pnt.hxx>
10 
11 #include <deal.II/base/point.h>
12 #include <deal.II/grid/tria.h>
13 #include <deal.II/fe/mapping.h>
15 #include <deal.II/lac/vector.h>
16 
17 using namespace dealii;
18 
19 #include "occ_normal_projection.h"
20 
21 namespace OpenCascade
22 {
23 
25  {
26  bool operator() (std::pair<double, double> a,std::pair<double, double> b) const
27  {
28  return a.second<b.second;
29  }
30  };
31 
32 
33 
35  {
36  public:
38  // euler_vector for which
39  // smoothing_dofs is set to
40  // true. The euler_vector is
41  // interpreted as a collection
42  // of triples associated with
43  // the displacement of the
44  // nodes of the mesh, like the
45  // one which is used in
46  // MappingQEulerian.
47  //
48  // A static Q1 mapping is used
49  // to compute the location of
50  // the reference
51  // points. Exceptions are
52  // thrown if the points are not
53  // located on the Geometry
54  // curve. The parameter
55  // base_point_id, identifies
56  // the 3 dofs associated with
57  // the base point, which is the
58  // one left untouched. All
59  // other points are moved
60  // proportionally to their
61  // arclength distance from the
62  // base point. The driving
63  // point id is given by
64  // drivint_point_id. These two
65  // ids are such that
66  // 3*base_point_id+i is the ith
67  // component of the base_point,
68  // contained in
69  // euler_vector(3*base_point_id+i).
70 
71  LineSmoothing(Vector<double> &euler_vector,
72  Handle(Geom_Curve) ref_curve,
73  TopLoc_Location *curr_loc,
74  const DoFHandler<2,3> &dh,
75  const std::vector<bool> &smoothing_dofs,
76  const unsigned int base_point_id,
77  const unsigned int driving_point_id,
78  const double tolerance = 1e-4);
79 
80  void update_reference(unsigned int base_point_id,
81  unsigned int driving_point_id);
96  void smooth(bool maintain_on_original_curve=true);
97 
98  inline Handle(Geom_Curve) get_curve()
99  {
100  return this->curve;
101  }
102 
104  {
105  return this->lengths_before_smoothing;
106  }
107 
109  {
110  return this->lengths_after_smoothing;
111  }
112 
113  inline std::vector<unsigned int> &get_node_indices()
114  {
115  return this->node_indices;
116  }
117 
118 
119  void get_curve_tangent_vectors_at_smoothing_dofs(Vector<double> &tangents);
120 
121  void get_curve_length_ratios_at_smoothing_dofs(Vector<double> &length_ratios);
122 
123 
124  typedef std::map<std::pair<double, double>, unsigned int, comp_points_on_curve >::iterator iterator;
125 
126  private:
128  Handle(Geom_Curve) ref_curve;
129  TopLoc_Location ref_location;
130  TopLoc_Location used_location;
131  TopLoc_Location *curr_location;
132  Handle(Geom_Curve) curve;
133  NormalProjection<1> projection;
134 
135  const DoFHandler<2,3> &dh;
136  const std::vector<bool> &smoothing_dofs;
137  unsigned int base_point_id;
138  unsigned int driving_point_id;
139  double occ_base_t;
140  double occ_driving_t;
141 
142  std::vector<Point<3> > support_points;
143 
144  double ref_L;
145  double tolerance;
146 
147  Vector<double> fixed_length_ratios;
148  Vector<double> lengths_before_smoothing;
149  Vector<double> lengths_after_smoothing;
150  std::vector<unsigned int> node_indices;
151  std::map<std::pair<double, double>, unsigned int, comp_points_on_curve > smoothing_list;
152  GeomLib_Tool tool;
153  };
154 
155 
156 }
157 
158 #endif
Vector< double > & get_lengths_after_smoothing()
std::vector< unsigned int > & get_node_indices()
Vector< double > & euler_vector
STL namespace.
Handle(Geom_Curve) get_curve()
std::map< std::pair< double, double >, unsigned int, comp_points_on_curve >::iterator iterator
We collect in this namespace all utilities which operate on OpenCascade entities which don't need cla...
Handle(Geom_Curve) NumericalTowingTank
Vector< double > & get_lengths_before_smoothing()