<h3>Specific improvements</h3>
-<ol>
+<ol>
+ <li> New: Add NURBSPatchManifold. This class is a child of ChartManifold and
+ implements a manifold descriptor for the face of a CAD imported usign
+ OpenCASCADE.
+ <br>
+ (Mauro Bardelloni, 2016/03/09)
+ </li>
+
<li> New: When using C++11, there is now a function Threads::new_task()
that can take as an argument either a lambda function, or the result
of a std::bind expression, or anything else that can be called as in a
};
/**
- * TODO:
+ * Manifold description for the face of a CAD imported usign OpenCASCADE.
*
* @ingroup manifold
*
- * @author Luca Heltai, 2014
+ * @author Andrea Mola, Mauro Bardelloni, 2016
*/
template <int dim, int spacedim>
class NURBSPatchManifold : public ChartManifold<dim, spacedim, 2>
{
public:
/**
- * TODO:
- * Asser: count shape!
+ * The constructor takes an OpenCASCADE TopoDS_Face @p face and an optional
+ * @p tolerance. This class uses the interval OpenCASCADE variables @var u,
+ * @var v to descrive the manifold.
*/
NURBSPatchManifold(const TopoDS_Face &face, const double tolerance = 1e-7);
/**
- * TODO: 3 -> 2
+ * 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;
/**
- * TODO:
+ * 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;
/**
* Given a point in the spacedim dimensional Euclidean space, this
* method returns the derivatives of the function $F$ that maps from
- * the polar coordinate system to the Euclidean coordinate
+ * the uv coordinate system to the Euclidean coordinate
* system. In other words, it is a matrix of size
- * $\text{spacedim}\times\text{spacedim}$.
+ * $\text{spacedim}\times\text{chartdim}$.
*
* This function is used in the computations required by the
* get_tangent_vector() function.
* Refer to the general documentation of this class for more information.
*/
virtual
- DerivativeForm<1,dim,spacedim>
+ DerivativeForm<1,2,spacedim>
push_forward_gradient(const Point<2> &chart_point) const;
- // /**
- // * Let the new point be the average sum of surrounding vertices.
- // *
- // * In the two dimensional implementation, we use the pull_back and
- // * push_forward mechanism. For three dimensions, this does not work well, so
- // * we overload the get_new_point function directly.
- // */
- // virtual Point<spacedim>
- // get_new_point(const Quadrature<spacedim> &quad) const;
+ private:
+ /**
+ * Return a tuple representing the minimum and maximum values of u
+ * and v. Precisely, it returns (u_min, u_max, v_min, v_max)
+ */
+ std_cxx11::tuple<double, double, double, double>
+ get_uv_bounds() const;
/**
- * TODO
+ * An OpenCASCADE TopoDS_Face @p face given by the CAD.
*/
- std::tuple<double, double, double, double>
- get_uv_bounds() const;
-
- private:
TopoDS_Face face;
+ /**
+ * Tolerance used by OpenCASCADE to identify points in each
+ * operation.
+ */
double tolerance;
};
surf->D1(chart_point[0],chart_point[1], q, Du, Dv);
DX[0][0] = Du.X();
- DX[0][1] = Du.Y();
- DX[0][2] = Du.Z();
- DX[1][0] = Dv.X();
+ DX[1][0] = Du.Y();
+ DX[2][0] = Du.Z();
+ DX[0][1] = Dv.X();
DX[1][1] = Dv.Y();
- DX[1][2] = Dv.Z();
+ DX[2][1] = Dv.Z();
return DX;
}
template<>
- std::tuple<double, double, double, double>
+ std_cxx11::tuple<double, double, double, double>
NURBSPatchManifold<2, 3>::
get_uv_bounds() const
{
Standard_Real umin, umax, vmin, vmax;
- BRepTools::UVBounds(face, umin, umax, vmin, vmax);
- return std::make_tuple(umin, umax, vmin, vmax);
+ BRepTools::UVBounds(face, umin, umax, vmin, vmax);
+ return std_cxx11::make_tuple(umin, umax, vmin, vmax);
}
-template class NURBSPatchManifold<2, 3 >;
// Explicit instantiations
+ template class NURBSPatchManifold<2, 3 >;
#include "boundary_lib.inst"
} // end namespace OpenCASCADE
dealii::OpenCASCADE::NURBSPatchManifold<2,3> manifold(faces[0]);
- auto bounds = manifold.get_uv_bounds();
- auto u_min = std::get<0>(bounds);
- auto u_max = std::get<1>(bounds);
- auto v_min = std::get<2>(bounds);
- auto v_max = std::get<3>(bounds);
-
- deallog << " u_min = " << std::get<0>(bounds) << std::endl;
- deallog << " u_max = " << std::get<1>(bounds) << std::endl;
- deallog << " v_min = " << std::get<2>(bounds) << std::endl;
- deallog << " v_max = " << std::get<3>(bounds) << std::endl;
+ const double u_min = 1.96115;
+ const double u_max = 10.0000;
+ const double v_min = 0.00000;
+ const double v_max = 11.0000;
deallog << "=======================================" << std::endl;
int len = 10;
- for(unsigned int i = 0; i <= 10; ++i)
- {
- double step = ((double)i)/10;
- deallog << " pos = " << step << std::endl;
- double u_pos = step * u_min + (1-step) * u_max;
- double v_pos = step * v_min + (1-step) * v_max;
- Point<2> uv(u_pos, v_pos);
- deallog << " uv = " << uv << std::endl;
- auto q = manifold.push_forward(uv);
- deallog << " q = " << q << std::endl;
- auto uv_ = manifold.pull_back(q);
- deallog << " uv = " << uv << std::endl;
- deallog << "=======================================" << std::endl;
- }
+ for (unsigned int i = 0; i <= 10; ++i)
+ {
+ double step = ((double)i)/10;
+ deallog << " pos = " << step << std::endl;
+ double u_pos = step * u_min + (1-step) * u_max;
+ double v_pos = step * v_min + (1-step) * v_max;
+ Point<2> uv(u_pos, v_pos);
+ deallog << " uv = " << uv << std::endl;
+ Point<3> q = manifold.push_forward(uv);
+ deallog << " q = " << q << std::endl;
+ Point<2> uv_ = manifold.pull_back(q);
+ deallog << " uv = " << uv << std::endl;
+ deallog << "=======================================" << std::endl;
+ }
double u_pos = (u_min + u_max)/2;
double v_pos = (v_min + v_max)/2;
Point<2> uv(u_pos, v_pos);
- auto D = manifold.push_forward_gradient(uv);
+ DerivativeForm<1,2,3> D = manifold.push_forward_gradient(uv);
deallog << "=======================================" << std::endl;
- deallog << " | " << D[1][0] << " | " << D[0][0] <<" | " << std::endl;
- deallog << " | " << D[1][1] << " | " << D[0][1] <<" | " << std::endl;
- deallog << " | " << D[1][2] << " | " << D[0][2] <<" | " << std::endl;
+ deallog << " | " << D[0][1] << " | " << D[0][0] <<" | " << std::endl;
+ deallog << " | " << D[1][1] << " | " << D[1][0] <<" | " << std::endl;
+ deallog << " | " << D[2][1] << " | " << D[2][0] <<" | " << std::endl;
deallog << "=======================================" << std::endl;
return 0;
-DEAL:: u_min = 1.96115
-DEAL:: u_max = 10.0000
-DEAL:: v_min = 0.00000
-DEAL:: v_max = 11.0000
DEAL::=======================================
DEAL:: pos = 0.00000
DEAL:: uv = 10.0000 11.0000
DEAL:: uv = 10.0000 11.0000
DEAL::=======================================
DEAL:: pos = 0.100000
-DEAL:: uv = 9.19611 9.90000
+DEAL:: uv = 9.19612 9.90000
DEAL:: q = -2.76096 0.0156017 0.279445
-DEAL:: uv = 9.19611 9.90000
+DEAL:: uv = 9.19612 9.90000
DEAL::=======================================
DEAL:: pos = 0.200000
DEAL:: uv = 8.39223 8.80000
DEAL:: uv = 8.39223 8.80000
DEAL::=======================================
DEAL:: pos = 0.300000
-DEAL:: uv = 7.58834 7.70000
+DEAL:: uv = 7.58835 7.70000
DEAL:: q = -2.53098 0.0178089 0.0300784
-DEAL:: uv = 7.58834 7.70000
+DEAL:: uv = 7.58835 7.70000
DEAL::=======================================
DEAL:: pos = 0.400000
DEAL:: uv = 6.78446 6.60000
-DEAL:: q = -2.42303 0.0237883 -0.0495104
+DEAL:: q = -2.42303 0.0237882 -0.0495104
DEAL:: uv = 6.78446 6.60000
DEAL::=======================================
DEAL:: pos = 0.500000
DEAL::=======================================
DEAL:: pos = 0.600000
DEAL:: uv = 5.17669 4.40000
-DEAL:: q = -2.21983 0.0279539 -0.128575
+DEAL:: q = -2.21983 0.0279540 -0.128575
DEAL:: uv = 5.17669 4.40000
DEAL::=======================================
DEAL:: pos = 0.700000
-DEAL:: uv = 4.37280 3.30000
+DEAL:: uv = 4.37281 3.30000
DEAL:: q = -2.12190 0.0130014 -0.163149
-DEAL:: uv = 4.37280 3.30000
+DEAL:: uv = 4.37281 3.30000
DEAL::=======================================
DEAL:: pos = 0.800000
DEAL:: uv = 3.56892 2.20000
-DEAL:: q = -2.04822 0.00230331 -0.181689
+DEAL:: q = -2.04822 0.00230333 -0.181689
DEAL:: uv = 3.56892 2.20000
DEAL::=======================================
DEAL:: pos = 0.900000
DEAL:: uv = 2.76503 1.10000
-DEAL:: q = -2.02052 0.00251434 -0.158299
+DEAL:: q = -2.02052 0.00251432 -0.158299
DEAL:: uv = 2.76503 1.10000
DEAL::=======================================
DEAL:: pos = 1.00000
DEAL:: uv = 1.96115 0.00000
-DEAL:: q = -1.98715 1.02653e-19 -0.149998
+DEAL:: q = -1.98715 1.02654e-19 -0.149998
DEAL:: uv = 1.96115 0.00000
DEAL::=======================================
DEAL::=======================================
-DEAL:: | -0.00164538 | -0.122937 |
-DEAL:: | -0.00311282 | -0.000314830 |
-DEAL:: | 0.0340939 | -0.00164538 |
+DEAL:: | -0.00164538 | -0.122938 |
+DEAL:: | -0.00311285 | -0.000314827 |
+DEAL:: | 0.0340939 | 8.33097e-05 |
DEAL::=======================================