diff --git a/dumux/multidomain/embedded/circlepoints.hh b/dumux/multidomain/embedded/circlepoints.hh
index 9c991d21a35d2fef786c82a9b5b5e71e9c411392..004115d391efcb26d4024a5bf6159ebdfbc14d53 100644
--- a/dumux/multidomain/embedded/circlepoints.hh
+++ b/dumux/multidomain/embedded/circlepoints.hh
@@ -82,21 +82,12 @@ void circlePoints(std::vector<GlobalPosition>& points,
 
 /*!
  * \ingroup EmbeddedCoupling
- * \brief returns a vector of points on a circle
- * \param center the circle center
- * \param normal the normal to the circle plane
- * \param radius the circle radius
+ * \brief returns a vector of sin and cos values of a circle parametrization
  * \param numPoints the number of points
  */
-template<class GlobalPosition, class Scalar>
-std::vector<GlobalPosition> circlePoints(const GlobalPosition& center,
-                                         const GlobalPosition& normal,
-                                         const Scalar radius,
-                                         const std::size_t numPoints = 20)
+template<class Scalar = double>
+std::vector<Scalar> circlePointsSinCos(const std::size_t numPoints)
 {
-    std::vector<GlobalPosition> points;
-
-    // precompute the sin/cos
     std::vector<Scalar> sincos(2*numPoints);
     Scalar t = 0 + 0.1; // add arbitrary offset
     for (std::size_t i = 0; i < numPoints; ++i)
@@ -107,7 +98,27 @@ std::vector<GlobalPosition> circlePoints(const GlobalPosition& center,
         t += 2*M_PI/numPoints;
         if(t > 2*M_PI) t -= 2*M_PI;
     }
+    return sincos;
+}
+
+/*!
+ * \ingroup EmbeddedCoupling
+ * \brief returns a vector of points on a circle
+ * \param center the circle center
+ * \param normal the normal to the circle plane
+ * \param radius the circle radius
+ * \param numPoints the number of points
+ */
+template<class GlobalPosition, class Scalar>
+std::vector<GlobalPosition> circlePoints(const GlobalPosition& center,
+                                         const GlobalPosition& normal,
+                                         const Scalar radius,
+                                         const std::size_t numPoints = 20)
+{
+    // precompute the sin/cos values
+    const auto sincos = circlePointsSinCos<Scalar>(numPoints);
 
+    std::vector<GlobalPosition> points;
     circlePoints(points, sincos, center, normal, radius);
     return points;
 }