NormalProjectionBoundary(const TopoDS_Shape &sh,
const double tolerance=1e-7);
+ /**
+ * Clone the current Manifold.
+ */
+ virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
/**
* Perform the actual projection onto the manifold. This function, in
* debug mode, checks that each of the @p surrounding_points is within
*/
virtual Point<spacedim>
project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
- const Point<spacedim> &candidate) const;
+ const Point<spacedim> &candidate) const override;
private:
const Tensor<1,spacedim> &direction,
const double tolerance=1e-7);
+ /**
+ * Clone the current Manifold.
+ */
+ virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
/**
* Perform the actual projection onto the manifold. This function, in
* debug mode, checks that each of the @p surrounding_points is within
*/
virtual Point<spacedim>
project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
- const Point<spacedim> &candidate) const;
+ const Point<spacedim> &candidate) const override;
private:
/**
NormalToMeshProjectionBoundary(const TopoDS_Shape &sh,
const double tolerance=1e-7);
+ /**
+ * Clone the current Manifold.
+ */
+ virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
/**
* Perform the actual projection onto the manifold. This function, in
* debug mode, checks that each of the @p surrounding_points is within
*/
virtual Point<spacedim>
project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
- const Point<spacedim> &candidate) const;
+ const Point<spacedim> &candidate) const override;
private:
/**
ArclengthProjectionLineManifold(const TopoDS_Shape &sh,
const double tolerance=1e-7);
+ /**
+ * Clone the current Manifold.
+ */
+ virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
/**
* Given a point on real space, find its arclength parameter. Throws an
* error in debug mode, if the point is not on the TopoDS_Edge given at
* construction time.
*/
virtual Point<1>
- pull_back(const Point<spacedim> &space_point) const;
+ pull_back(const Point<spacedim> &space_point) const override;
/**
* Given an arclength parameter, find its image in real space.
*/
virtual Point<spacedim>
- push_forward(const Point<1> &chart_point) const;
+ push_forward(const Point<1> &chart_point) const override;
private:
+ /**
+ * The actual shape used to build this object.
+ */
+ const TopoDS_Shape sh;
+
/**
* A Curve adaptor. This is the one which is used in the computations, and
* it points to the right one above.
*/
NURBSPatchManifold(const TopoDS_Face &face, const double tolerance = 1e-7);
+ /**
+ * Clone the current Manifold.
+ */
+ virtual std::unique_ptr<Manifold<dim,spacedim> > clone() const override;
+
/**
* Pull back the given point from the Euclidean space. Will return the uv
* coordinates associated with the point @p space_point.
*/
virtual Point<2>
- pull_back(const Point<spacedim> &space_point) const;
+ pull_back(const Point<spacedim> &space_point) const override;
/**
* Given a @p chart_point in the uv coordinate system, this method returns the
* Euclidean coordinates associated.
*/
virtual Point<spacedim>
- push_forward(const Point<2> &chart_point) const;
+ push_forward(const Point<2> &chart_point) const override;
/**
* Given a point in the spacedim dimensional Euclidean space, this
*/
virtual
DerivativeForm<1,2,spacedim>
- push_forward_gradient(const Point<2> &chart_point) const;
+ push_forward_gradient(const Point<2> &chart_point) const override;
private:
/**
}
+
+ template<int dim, int spacedim>
+ std::unique_ptr<Manifold<dim, spacedim> >
+ NormalProjectionBoundary<dim,spacedim>::clone() const
+ {
+ return std::unique_ptr<Manifold<dim,spacedim> >(new NormalProjectionBoundary(sh, tolerance));
+ }
+
+
+
template <int dim, int spacedim>
Point<spacedim> NormalProjectionBoundary<dim,spacedim>::
project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
}
+
+ template<int dim, int spacedim>
+ std::unique_ptr<Manifold<dim, spacedim> >
+ DirectionalProjectionBoundary<dim,spacedim>::clone() const
+ {
+ return std::unique_ptr<Manifold<dim,spacedim> >
+ (new DirectionalProjectionBoundary(sh, direction, tolerance));
+ }
+
+
+
template <int dim, int spacedim>
Point<spacedim> DirectionalProjectionBoundary<dim,spacedim>::
project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
ExcMessage("NormalToMeshProjectionBoundary needs a shape containing faces to operate."));
}
+ template<int dim, int spacedim>
+ std::unique_ptr<Manifold<dim, spacedim> >
+ NormalToMeshProjectionBoundary<dim,spacedim>::clone() const
+ {
+ return std::unique_ptr<Manifold<dim, spacedim> >
+ (new NormalToMeshProjectionBoundary<dim,spacedim>(sh,tolerance));
+ }
+
template <int dim, int spacedim>
Point<spacedim> NormalToMeshProjectionBoundary<dim,spacedim>::
ChartManifold<dim,spacedim,1>(sh.Closed() ?
Point<1>(shape_length(sh)) :
Point<1>()),
+ sh(sh),
curve(curve_adaptor(sh)),
tolerance(tolerance),
length(shape_length(sh))
}
+
+ template<int dim, int spacedim>
+ std::unique_ptr<Manifold<dim, spacedim> >
+ ArclengthProjectionLineManifold<dim, spacedim>::clone() const
+ {
+ return std::unique_ptr<Manifold<dim, spacedim> >
+ (new ArclengthProjectionLineManifold(sh,tolerance));
+ }
+
+
+
template <int dim, int spacedim>
Point<1>
ArclengthProjectionLineManifold<dim,spacedim>::pull_back(const Point<spacedim> &space_point) const
tolerance(tolerance)
{}
+
+
+ template<int dim, int spacedim>
+ std::unique_ptr<Manifold<dim, spacedim> >
+ NURBSPatchManifold<dim,spacedim>::clone() const
+ {
+ return std::unique_ptr<Manifold<dim, spacedim> >
+ (new NURBSPatchManifold<dim,spacedim>(face,tolerance));
+ }
+
+
+
template <int dim, int spacedim> Point<2>
NURBSPatchManifold<dim, spacedim>::
pull_back(const Point<spacedim> &space_point) const