From: Rene Gassmoeller Date: Sun, 3 Dec 2017 21:19:42 +0000 (-0700) Subject: Address some comments X-Git-Tag: v9.0.0-rc1~626^2~5 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=9a6a7be5249e5d15e5163ebdc44b1311e685d39c;p=dealii.git Address some comments --- diff --git a/include/deal.II/grid/manifold_lib.h b/include/deal.II/grid/manifold_lib.h index 59c38872a1..980323d031 100644 --- a/include/deal.II/grid/manifold_lib.h +++ b/include/deal.II/grid/manifold_lib.h @@ -253,7 +253,7 @@ public: void get_new_points (const ArrayView> &surrounding_points, const Table<2,double> &weights, - ArrayView> new_points) const; + ArrayView> new_points) const override; /** * Return a point on the spherical manifold which is intermediate @@ -307,8 +307,6 @@ private: const PolarManifold polar_manifold; }; - - /** * Cylindrical Manifold description. In three dimensions, points are * transformed using a cylindrical coordinate system along the x-, diff --git a/source/grid/manifold_lib.cc b/source/grid/manifold_lib.cc index c3535ab51d..b66bb88043 100644 --- a/source/grid/manifold_lib.cc +++ b/source/grid/manifold_lib.cc @@ -413,7 +413,7 @@ get_new_points (const ArrayView> &surrounding_points, // Step 1: Check for some special cases, create simple linear guesses otherwise. const double tolerance = 1e-10; - std::vector accurate_point_was_found(new_points.size(),false); + boost::container::small_vector accurate_point_was_found(new_points.size(),false); const ArrayView> array_directions = make_array_view(directions.begin(),directions.end()); const ArrayView array_distances = make_array_view(distances.begin(),distances.end()); for (unsigned int row=0; row> &surrounding_points, // Search for duplicate directions and merge them to minimize the cost of // the get_new_point function call below. - Table<2,double> merged_weights(weights); + boost::container::small_vector merged_weights(weights.size(0)*weights.size(1)); boost::container::small_vector, 100> merged_directions(surrounding_points.size(),Point()); boost::container::small_vector merged_distances(surrounding_points.size(),0.0); @@ -464,7 +464,7 @@ get_new_points (const ArrayView> &surrounding_points, { bool found_duplicate = false; - // This inner loop if of N^2 complexity, but surrounding_points.size() + // This inner loop is of $O(N^2)$ complexity, but surrounding_points.size() // is usually at most 8 points large. for (unsigned int j = 0; j < n_unique_directions; ++j) { @@ -473,7 +473,7 @@ get_new_points (const ArrayView> &surrounding_points, { found_duplicate = true; for (unsigned int row = 0; row < weights.size(0); ++row) - merged_weights[row][j] += weights[row][i]; + merged_weights[row*weights.size(1) + j] += weights[row][i]; } } @@ -482,7 +482,7 @@ get_new_points (const ArrayView> &surrounding_points, merged_directions[n_unique_directions] = directions[i]; merged_distances[n_unique_directions] = distances[i]; for (unsigned int row = 0; row < weights.size(0); ++row) - merged_weights[row][n_unique_directions] = weights[row][i]; + merged_weights[row*weights.size(1) + n_unique_directions] = weights[row][i]; ++n_unique_directions; } @@ -497,7 +497,7 @@ get_new_points (const ArrayView> &surrounding_points, for (unsigned int row=0; row array_merged_weights (merged_weights[row].begin(),n_unique_directions); + const ArrayView array_merged_weights (&merged_weights[row*weights.size(1)],n_unique_directions); new_points[row] = get_new_point(array_merged_directions, array_merged_distances, array_merged_weights, @@ -560,6 +560,141 @@ guess_new_point(const ArrayView> &directions, } +namespace +{ + template + Point + do_get_new_point(const ArrayView> &directions, + const ArrayView &distances, + const ArrayView &weights, + const Point candidate_point, + const Point center) + { + Assert(false,ExcNotImplemented()); + } + + template <> + Point<3> + do_get_new_point(const ArrayView> &directions, + const ArrayView &distances, + const ArrayView &weights, + const Point<3> candidate_point, + const Point<3> center) + { + AssertDimension(directions.size(), distances.size()); + AssertDimension(directions.size(), weights.size()); + + const unsigned int n_merged_points = directions.size(); + const double tolerance = 1e-10; + const int max_iterations = 10; + + // Recover radius and normalized direction from candidate point + Tensor <1,3> candidate = candidate_point - center; + const double rho = candidate.norm(); + candidate /= rho; + + { + // If the candidate happens to coincide with a normalized + // direction, we return it. Otherwise, the Hessian would be singular. + for (unsigned int i=0; i unit_manifold; + Assert(std::abs(weights[0] + weights[1] - 1.0) < 1e-13, + ExcMessage("Weights do not sum up to 1")); + Point<3> intermediate = unit_manifold.get_intermediate_point + (Point<3>(directions[0]), Point<3>(directions[1]), weights[1]); + return center + rho * intermediate; + } + + Tensor<1,3> vPerp; + Tensor<2,2> Hessian; + Tensor<1,2> gradient; + Tensor<1,2> gradlocal; + + // On success we exit the loop early. + // Otherwise, we just take the result after max_iterations steps. + for (unsigned int i=0; i Clocalx = internal::compute_normal(candidate); + const Tensor<1,3> Clocaly = cross_product_3d(candidate, Clocalx); + + // For each vertices vector, compute the tangent vector from candidate + // towards the vertices vector -- its length is the spherical length + // from candidate to the vertices vector. + // Then compute its contribution to the Hessian. + gradient = 0.; + Hessian = 0.; + for (unsigned int i=0; i1.e-15) + { + vPerp = internal::projected_direction(directions[i], candidate); + const double sintheta = vPerp.norm(); + if (sinthetatolerance, ExcInternalError()); + + const Tensor<2,2> inverse_Hessian = invert(Hessian); + + const Tensor<1,2> xDisplocal = inverse_Hessian*gradient; + const Tensor<1,3> xDisp = xDisplocal[0]*Clocalx + xDisplocal[1]*Clocaly; + + // Step 2b: rotate candidate in direction xDisp for a new candidate. + const Tensor<1,3> candidateOld = candidate; + candidate = internal::apply_exponential_map(candidate, xDisp); + + // Step 3c: return the new candidate if we didn't move + if ((candidate-candidateOld).norm_square() < tolerance*tolerance) + break; + } + + Assert (std::abs(candidate[2]) < tolerance, + ExcInternalError()); + } + return center + rho*candidate; + } +} + + template Point @@ -569,142 +704,47 @@ get_new_point (const ArrayView> &directions, const ArrayView &weights, const Point candidate_point) const { - AssertDimension(directions.size(), distances.size()); - AssertDimension(directions.size(), weights.size()); - - const unsigned int n_merged_points = directions.size(); - const double tolerance = 1e-10; - const int max_iterations = 10; - - // Recover radius and normalized direction from candidate point - Tensor <1,spacedim> candidate = candidate_point - center; - const double rho = candidate.norm(); - candidate /= rho; - - // In this step, we consider all points and directions to be embedded - // in a three-dimensional space. The case spacedim < 2 was handled in get_new_points() - { - Tensor<1, 3> xVec; - for (unsigned int c=0; c, 100> directions_3d(directions.size()); - for (unsigned int i=0; i normalized_direction = directions[i] / distances[i]; - const double squared_distance = (candidate - normalized_direction).norm_square(); - if (squared_distance < tolerance*tolerance) - return center + rho * candidate; - - // append direction. check if the normalized candidate direction is - // the same as a previous direction (to a tighter tolerance (1e-14)^2 - // than the outer ones to really not miss anything) -> in that case we - // can simply add the weights. Since the trigonometric functions used - // below are quite expensive, it makes sense to merge the points here, - // even if this search loop is of quadratic complexity loop (but we - // rarely have more than 9 points) - Tensor<1,3> direction_3d; - for (unsigned int c=0; c unit_manifold; - Assert(std::abs(weights[0] + weights[1] - 1.0) < 1e-13, - ExcMessage("Weights do not sum up to 1")); - Point<3> intermediate = unit_manifold.get_intermediate_point - (Point<3>(directions_3d[0]), Point<3>(directions_3d[1]), weights[1]); - // copy back to spacedim-point - Point p; - for (unsigned int d=0; d(); +} - Tensor<1,3> vPerp; - Tensor<2,2> Hessian; - Tensor<1,2> gradient; - Tensor<1,2> gradlocal; - // On success we exit the loop early. - // Otherwise, we just take the result after max_iterations steps. - for (unsigned int i=0; i Clocalx = internal::compute_normal(xVec); - const Tensor<1,3> Clocaly = cross_product_3d(xVec, Clocalx); - - // For each vertices vector, compute the tangent vector from xVec - // towards the vertices vector -- its length is the spherical length - // from xVec to the vertices vector. - // Then compute its contribution to the Hessian. - gradient = 0.; - Hessian = 0.; - for (unsigned int i=0; i1.e-15) - { - vPerp = internal::projected_direction(directions_3d[i], xVec); - const double sintheta = vPerp.norm(); - if (sinthetatolerance, ExcInternalError()); +template <> +Point<3> +SphericalManifold<1,3>:: +get_new_point (const ArrayView> &directions, + const ArrayView &distances, + const ArrayView &weights, + const Point<3> candidate_point) const +{ + return do_get_new_point(directions,distances,weights,candidate_point,center); +} - const Tensor<2,2> inverse_Hessian = invert(Hessian); - const Tensor<1,2> xDisplocal = inverse_Hessian*gradient; - const Tensor<1,3> xDisp = xDisplocal[0]*Clocalx + xDisplocal[1]*Clocaly; - // Step 2b: rotate xVec in direction xDisp for a new candidate. - const Tensor<1,3> xVecOld = xVec; - xVec = internal::apply_exponential_map(xVec, xDisp); +template <> +Point<3> +SphericalManifold<2,3>:: +get_new_point (const ArrayView> &directions, + const ArrayView &distances, + const ArrayView &weights, + const Point<3> candidate_point) const +{ + return do_get_new_point(directions,distances,weights,candidate_point,center); +} - // Step 3c: return the new candidate if we didn't move - if ((xVec-xVecOld).norm_square() < tolerance*tolerance) - break; - } - for (unsigned int c=0; c +Point<3> +SphericalManifold<3,3>:: +get_new_point (const ArrayView> &directions, + const ArrayView &distances, + const ArrayView &weights, + const Point<3> candidate_point) const +{ + return do_get_new_point(directions,distances,weights,candidate_point,center); }