3#include <dune/common/fmatrix.hh>
10namespace ProjectionImplementation {
22template<
typename Coordinate,
typename Field>
26 Coordinate x(Field(0));
42inline std::pair<unsigned, unsigned>
43edgeToCorners(
unsigned edge)
46 case 0:
return {0, 1};
47 case 1:
return {0, 2};
48 case 2:
return {1, 2};
50 DUNE_THROW(Dune::Exception,
"Unexpected edge number.");
68template<
typename Coordinate,
typename Corners>
69inline typename Corners::value_type
70interpolate(
const Coordinate& x,
const Corners& corners)
73 for (
unsigned i = 0; i < corners.size() - 1; ++i)
74 y.axpy(x[i], corners[i+1] - corners[0]);
89template<
typename Coordinate,
typename Normals>
90inline typename Normals::value_type
91interpolate_unit_normals(
const Coordinate& x,
const Normals& normals)
93 auto n = interpolate(x, normals);
109template<
typename Coordinate,
typename Field>
111inside(
const Coordinate& x,
const Field& epsilon)
113 const unsigned dim = Coordinate::dimension;
115 for (
unsigned i = 0; i < dim-1; ++i) {
121 if (sum <= Field(1) + epsilon)
128template<
typename Coordinate>
132 , m_max_normal_product(max_normal_product)
137template<
typename Coordinate>
145template<
typename Coordinate>
146template<
typename Corners,
typename Normals>
149::doProjection(
const std::tuple<Corners&, Corners&>& corners,
const std::tuple<Normals&, Normals&>& normals)
162 using namespace ProjectionImplementation;
164 typedef Dune::FieldMatrix<Field, dim, dim> Matrix;
167 const auto& origin = get<0>(corners);
168 const auto& origin_normals = get<0>(normals);
169 const auto& target = get<1>(corners);
170 const auto& target_normals = get<1>(normals);
171 auto& images = get<0>(m_images);
172 auto& success = get<0>(m_success);
178 std::array<Coordinate, dim-1> directions;
179 std::array<Field, dim-1> scales;
182 for (
unsigned i = 0; i < dim-1; ++i) {
183 directions[i] = target[i+1] - target[0];
184 scales[i] = directions[i].infinity_norm();
185 directions[i] /= scales[i];
186 scaleSum += scales[i];
189 for (
unsigned i = 0; i < dim-1; ++i) {
190 for (
unsigned j = 0; j < dim; ++j) {
191 m[j][i] = directions[i][j];
195 m_projection_valid =
true;
199 for (
unsigned i = 0; i < origin.size(); ++i) {
200 for (
unsigned j = 0; j < dim; ++j)
201 m[j][dim-1] = origin_normals[i][j];
203 const Coordinate rhs = origin[i] - target[0];
209 for (
unsigned j = 0; j < dim-1; ++j)
212 y[dim-1] *= Field(-1);
221 if(y[dim-1] < -2*scaleSum) {
222 success.set(i,
false);
223 m_projection_valid =
false;
227 const bool feasible = projectionFeasible(origin[i], origin_normals[i], y, target, target_normals);
228 success.set(i, feasible);
230 catch (
const Dune::FMatrixError&) {
231 success.set(i,
false);
232 m_projection_valid =
false;
237template<
typename Coordinate>
238template<
typename Corners,
typename Normals>
240Projection<Coordinate>
241::doInverseProjection(
const std::tuple<Corners&, Corners&>& corners,
const std::tuple<Normals&, Normals&>& normals)
253 using namespace ProjectionImplementation;
255 typedef Dune::FieldMatrix<Field, dim-1, dim-1> Matrix;
256 typedef Dune::FieldVector<Field, dim-1> Vector;
261 if (!m_projection_valid) {
262 get<1>(m_success).reset();
266 const auto& images = get<0>(m_images);
267 const auto& target_corners = get<1>(corners);
268 auto& preimages = get<1>(m_images);
269 auto& success = get<1>(m_success);
271 std::array<Coordinate, dim> v;
272 for (
unsigned i = 0; i < dim-1; ++i) {
273 v[i] = interpolate(images[i+1], target_corners);
274 v[i] -= interpolate(images[0], target_corners);
278 for (
unsigned i = 0; i < dim-1; ++i) {
279 for (
unsigned j = 0; j < dim-1; ++j) {
284 for (
unsigned i = 0; i < dim; ++i) {
286 v[dim-1] = target_corners[i];
287 v[dim-1] -= interpolate(images[0], target_corners);
290 for (
unsigned j = 0; j < dim-1; ++j)
291 rhs[j] = v[dim-1]*v[j];
294 for (
unsigned j = 0; j < dim-1; ++j)
295 preimages[i][j] = z[j];
298 const auto x = interpolate(z, get<0>(corners));
299 preimages[i][dim-1] = (x - target_corners[i]) * get<1>(normals)[i];
302 const bool feasible = projectionFeasible(target_corners[i], get<1>(normals)[i], preimages[i], get<0>(corners), get<0>(normals));
303 success.set(i, feasible);
307template<
typename Coordinate>
308template<
typename Corners,
typename Normals>
310Projection<Coordinate>
311::doEdgeIntersection(
const std::tuple<Corners&, Corners&>& corners,
const std::tuple<Normals&, Normals&>& normals)
313 using namespace ProjectionImplementation;
316 m_number_of_edge_intersections = 0;
327 if (!m_projection_valid || get<0>(m_success).all() || get<1>(m_success).all()) {
331 const auto& images = get<0>(m_images);
332 const auto& ys = get<1>(corners);
343 for (
unsigned edgex = 0; edgex < dim; ++edgex) {
345 std::tie(i, j) = edgeToCorners(edgex);
348 if (get<0>(m_success)[i] && get<0>(m_success)[j])
351 const auto pxi = interpolate(images[i], ys);
352 const auto pxj = interpolate(images[j], ys);
353 const auto pxjpxi = pxj - pxi;
355 typedef Dune::FieldMatrix<Field, dim-1, dim-1> Matrix;
356 typedef Dune::FieldVector<Field, dim-1> Vector;
358 for (
unsigned edgey = 0; edgey < dim; ++edgey) {
360 std::tie(k, l) = edgeToCorners(edgey);
363 if (get<1>(m_success)[k] && get<1>(m_success)[l])
366 const auto ykyl = ys[k] - ys[l];
367 const auto ykpxi = ys[k] - pxi;
370 bool parallel =
true;
371 for (
unsigned h=0; h<3; h++)
372 parallel &= std::abs(ykyl[(h+1)%3]*pxjpxi[(h+2)%3] - ykyl[(h+2)%3]*pxjpxi[(h+1)%3])<1e-14;
379 for (
unsigned m = 0; m < dim-1; ++m) {
380 const auto ym1y0 = ys[m+1] - ys[0];
381 mat[m][0] = pxjpxi * ym1y0;
382 mat[m][1] = ykyl * ym1y0;
383 rhs[m] = ykpxi * ym1y0;
392 if (!isfinite(z[0]) || !isfinite(z[1]))
396 if (z[0] < m_epsilon || z[0] > Field(1) - m_epsilon
397 || z[1] < m_epsilon || z[1] > Field(1) - m_epsilon)
400 Coordinate local_x = corner<Coordinate, Field>(i);
401 local_x.axpy(z[0], corner<Coordinate, Field>(j) - corner<Coordinate, Field>(i));
402 Coordinate local_y = corner<Coordinate, Field>(k);
403 local_y.axpy(z[1], corner<Coordinate, Field>(l) - corner<Coordinate, Field>(k));
406 if (!inside(local_x, m_epsilon) || !inside(local_y, m_epsilon))
410 auto xy = interpolate(local_x, get<0>(corners));
411 xy -= interpolate(local_y, get<1>(corners));
412 const auto nx = interpolate_unit_normals(local_x, get<0>(normals));
413 const auto ny = interpolate_unit_normals(local_y, get<1>(normals));
414 local_x[dim-1] = -(xy*nx);
415 local_y[dim-1] = xy*ny;
417 if (local_x[dim-1] < -m_overlap-m_epsilon || local_y[dim-1] < -m_overlap-m_epsilon)
421 if (nx*ny > m_max_normal_product + m_epsilon)
425 auto& intersection = m_edge_intersections[m_number_of_edge_intersections++];
426 intersection = { {{edgex, edgey}}, {{local_x, local_y}} };
428 catch(
const Dune::FMatrixError&) {
435template<
typename Coordinate>
436template<
typename Corners,
typename Normals>
437bool Projection<Coordinate>
438::projectionFeasible(
const Coordinate& x,
const Coordinate& nx,
const Coordinate& px,
const Corners& corners,
const Normals& normals)
const
440 using namespace ProjectionImplementation;
443 if (!inside(px, m_epsilon))
447 if (px[dim-1] < -m_overlap-m_epsilon)
452 xmy -= interpolate(px, corners);
453 const auto n = interpolate_unit_normals(px, normals);
454 const auto d = xmy * n;
455 if (d < -m_overlap-m_epsilon)
459 if (nx * n > m_max_normal_product + m_epsilon)
466template<
typename Coordinate>
467template<
typename Corners,
typename Normals>
469::project(
const std::tuple<Corners&, Corners&>& corners,
const std::tuple<Normals&, Normals&>& normals)
471 doProjection(corners, normals);
472 doInverseProjection(corners, normals);
473 doEdgeIntersection(corners, normals);
Projection of a line (triangle) on another line (triangle).
Definition: projection.hh:21
Coordinate::field_type Field
Scalar type.
Definition: projection.hh:61
Projection(const Field overlap=Field(0), const Field max_normal_product=Field(-0.1))
Definition: projection_impl.hh:130
void epsilon(const Field epsilon)
Set epsilon used for floating-point comparisons.
Definition: projection_impl.hh:140
void project(const std::tuple< Corners &, Corners & > &corners, const std::tuple< Normals &, Normals & > &normals)
Do the actual projection.
Definition: projection_impl.hh:469