+ /**
+ * Compute the arc cosine of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{acos(x[0]), acos(x[1]), ...,
+ * acos(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ acos(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::acos(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the arc sine of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{asin(x[0]), asin(x[1]), ...,
+ * asin(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ asin(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::asin(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the arc tangent of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{atan(x[0]), atan(x[1]), ...,
+ * atan(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ atan(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::atan(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the hyperbolic cosine of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{cosh(x[0]), cosh(x[1]), ...,
+ * cosh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ cosh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::cosh(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the hyperbolic sine of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{sinh(x[0]), sinh(x[1]), ...,
+ * sinh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ sinh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::sinh(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the hyperbolic tangent of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{tanh(x[0]), tanh(x[1]), ...,
+ * tanh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ tanh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::tanh(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the area hyperbolic cosine of a vectorized data field. The result
+ * is returned as vectorized array in the form <tt>{acosh(x[0]), acosh(x[1]),
+ * ..., acosh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ acosh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::acosh(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the area hyperbolic sine of a vectorized data field. The result is
+ * returned as vectorized array in the form <tt>{asinh(x[0]), asinh(x[1]),
+ * ..., asinh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ asinh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::asinh(x[i]);
+ return out;
+ }
+
+
+
+ /**
+ * Compute the area hyperbolic tangent of a vectorized data field. The result
+ * is returned as vectorized array in the form <tt>{atanh(x[0]), atanh(x[1]),
+ * ..., atanh(x[size()-1])}</tt>.
+ *
+ * @relatesalso VectorizedArray
+ */
+ template <typename Number, std::size_t width>
+ inline ::dealii::VectorizedArray<Number, width>
+ atanh(const ::dealii::VectorizedArray<Number, width> &x)
+ {
+ ::dealii::VectorizedArray<Number, width> out;
+ for (unsigned int i = 0; i < dealii::VectorizedArray<Number, width>::size();
+ ++i)
+ out[i] = std::atanh(x[i]);
+ return out;
+ }
+
+
+
/**
* Compute the exponential of a vectorized data field. The result is
* returned as vectorized array in the form <tt>{exp(x[0]), exp(x[1]), ...,
AssertThrow(std::fabs(e[i] - std::sin(d[i])) <
10. * std::numeric_limits<Number>::epsilon(),
ExcInternalError());
+ deallog << "OK" << std::endl << "Arc sine: ";
+ d = std::asin(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(d[i] - std::asin(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
deallog << "OK" << std::endl << "Cosine: ";
e = std::cos(c);
for (unsigned int i = 0; i < n_vectors; ++i)
AssertThrow(std::fabs(e[i] - std::cos(c[i])) <
10. * std::numeric_limits<Number>::epsilon(),
ExcInternalError());
+ deallog << "OK" << std::endl << "Arc cosine: ";
+ c = std::acos(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(c[i] - std::acos(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
deallog << "OK" << std::endl << "Tangent: ";
d = std::tan(e);
for (unsigned int i = 0; i < n_vectors; ++i)
AssertThrow(std::fabs(d[i] - std::tan(e[i])) <
10. * std::numeric_limits<Number>::epsilon(),
ExcInternalError());
+ deallog << "OK" << std::endl << "Arc tangent: ";
+ d = std::atan(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(d[i] - std::atan(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Hyperbolic cosine: ";
+ e = std::cosh(c);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(e[i] - std::cosh(c[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Area hyperbolic cosine: ";
+ c = std::acosh(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(c[i] - std::acosh(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Hyperbolic sine: ";
+ e = std::sinh(d);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(e[i] - std::sinh(d[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Area hyperbolic sine: ";
+ d = std::asinh(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(d[i] - std::asinh(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Hyperbolic tangent: ";
+ d = std::tanh(e);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(d[i] - std::tanh(e[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
+ deallog << "OK" << std::endl << "Area hyperbolic tangent: ";
+ e = std::atanh(d);
+ for (unsigned int i = 0; i < n_vectors; ++i)
+ AssertThrow(std::fabs(e[i] - std::atanh(d[i])) <
+ 10. * std::numeric_limits<Number>::epsilon(),
+ ExcInternalError());
deallog << "OK" << std::endl << "Exponential: ";
d = std::exp(c - a);
for (unsigned int i = 0; i < n_vectors; ++i)
DEAL:double::Maximum value: OK
DEAL:double::Minimum value: OK
DEAL:double::Sine: OK
+DEAL:double::Arc sine: OK
DEAL:double::Cosine: OK
+DEAL:double::Arc cosine: OK
DEAL:double::Tangent: OK
+DEAL:double::Arc tangent: OK
+DEAL:double::Hyperbolic cosine: OK
+DEAL:double::Area hyperbolic cosine: OK
+DEAL:double::Hyperbolic sine: OK
+DEAL:double::Area hyperbolic sine: OK
+DEAL:double::Hyperbolic tangent: OK
+DEAL:double::Area hyperbolic tangent: OK
DEAL:double::Exponential: OK
DEAL:double::Logarithm: OK
DEAL:float::Addition: OK
DEAL:float::Maximum value: OK
DEAL:float::Minimum value: OK
DEAL:float::Sine: OK
+DEAL:float::Arc sine: OK
DEAL:float::Cosine: OK
+DEAL:float::Arc cosine: OK
DEAL:float::Tangent: OK
+DEAL:float::Arc tangent: OK
+DEAL:float::Hyperbolic cosine: OK
+DEAL:float::Area hyperbolic cosine: OK
+DEAL:float::Hyperbolic sine: OK
+DEAL:float::Area hyperbolic sine: OK
+DEAL:float::Hyperbolic tangent: OK
+DEAL:float::Area hyperbolic tangent: OK
DEAL:float::Exponential: OK
DEAL:float::Logarithm: OK
DEAL:long double::Addition: OK
DEAL:long double::Maximum value: OK
DEAL:long double::Minimum value: OK
DEAL:long double::Sine: OK
+DEAL:long double::Arc sine: OK
DEAL:long double::Cosine: OK
+DEAL:long double::Arc cosine: OK
DEAL:long double::Tangent: OK
+DEAL:long double::Arc tangent: OK
+DEAL:long double::Hyperbolic cosine: OK
+DEAL:long double::Area hyperbolic cosine: OK
+DEAL:long double::Hyperbolic sine: OK
+DEAL:long double::Area hyperbolic sine: OK
+DEAL:long double::Hyperbolic tangent: OK
+DEAL:long double::Area hyperbolic tangent: OK
DEAL:long double::Exponential: OK
DEAL:long double::Logarithm: OK