20 #ifndef LIBMESH_TYPE_VECTOR_H 21 #define LIBMESH_TYPE_VECTOR_H 37 template <
typename T>
class TypeTensor;
38 template <
typename T>
class VectorValue;
39 template <
typename T>
class TensorValue;
52 template <
typename T2>
77 template <
typename Scalar1,
typename Scalar2,
typename Scalar3>
80 const Scalar1>::type x,
83 const Scalar2>::type y=0,
86 const Scalar3>::type z=0);
94 template <
typename Scalar>
98 const Scalar>::type * sfinae =
nullptr);
115 template <
typename T2>
126 template <
typename T2>
132 template <
typename Scalar>
137 { libmesh_assert_equal_to (p, Scalar(0)); this->
zero();
return *
this; }
142 const T &
operator () (
const unsigned int i)
const;
143 const T &
slice (
const unsigned int i)
const {
return (*
this)(i); }
149 T &
slice (
const unsigned int i) {
return (*
this)(i); }
156 template <
typename T2>
165 template <
typename T2>
171 template <
typename T2>
177 template <
typename T2>
185 template <
typename T2>
194 template <
typename T2>
200 template <
typename T2>
207 template <
typename T2>
220 template <
typename Scalar>
238 template <
typename Scalar>
257 template <
typename T2>
264 template <
typename T2>
271 template <
typename T2>
286 #ifdef LIBMESH_ENABLE_DEPRECATED 302 #ifdef LIBMESH_ENABLE_DEPRECATED 351 bool operator < (const TypeVector<T> & rhs)
const;
359 bool operator <= (const TypeVector<T> & rhs)
const;
390 friend std::ostream & operator << (std::ostream & os, const TypeVector<T> & t)
419 template <
typename T>
436 template <
typename T>
447 libmesh_assert_equal_to (y, 0);
453 libmesh_assert_equal_to (z, 0);
458 template <
typename T>
459 template <
typename Scalar1,
typename Scalar2,
typename Scalar3>
463 const Scalar1>::type x,
466 const Scalar2>::type y,
469 const Scalar3>::type z)
476 libmesh_assert_equal_to (y, 0);
482 libmesh_assert_equal_to (z, 0);
488 template <
typename T>
489 template <
typename Scalar>
494 const Scalar>::type * )
509 template <
typename T>
510 template <
typename T2>
515 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
516 _coords[i] = p._coords[i];
521 template <
typename T>
529 template <
typename T>
530 template <
typename T2>
534 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
535 _coords[i] = p._coords[i];
540 template <
typename T>
544 libmesh_assert_less (i, LIBMESH_DIM);
551 template <
typename T>
555 libmesh_assert_less (i, LIBMESH_DIM);
562 template <
typename T>
563 template <
typename T2>
575 _coords[1] + p._coords[1]);
580 _coords[1] + p._coords[1],
581 _coords[2] + p._coords[2]);
588 template <
typename T>
589 template <
typename T2>
600 template <
typename T>
601 template <
typename T2>
606 _coords[0] += p._coords[0];
610 _coords[0] += p._coords[0];
611 _coords[1] += p._coords[1];
615 _coords[0] += p._coords[0];
616 _coords[1] += p._coords[1];
617 _coords[2] += p._coords[2];
624 template <
typename T>
625 template <
typename T2>
630 _coords[0] += factor*p(0);
634 _coords[0] += factor*p(0);
635 _coords[1] += factor*p(1);
639 _coords[0] += factor*p(0);
640 _coords[1] += factor*p(1);
641 _coords[2] += factor*p(2);
648 template <
typename T>
649 template <
typename T2>
662 _coords[1] - p._coords[1]);
667 _coords[1] - p._coords[1],
668 _coords[2] - p._coords[2]);
675 template <
typename T>
676 template <
typename T2>
687 template <
typename T>
688 template <
typename T2>
692 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
693 _coords[i] -= p._coords[i];
698 template <
typename T>
699 template <
typename T2>
703 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
704 _coords[i] -= factor*p(i);
709 template <
typename T>
733 template <
typename T>
734 template <
typename Scalar>
761 template <
typename T,
typename Scalar>
774 template <
typename T>
779 _coords[0] *= factor;
783 _coords[0] *= factor;
784 _coords[1] *= factor;
788 _coords[0] *= factor;
789 _coords[1] *= factor;
790 _coords[2] *= factor;
798 template <
typename T>
799 template <
typename Scalar>
806 libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
830 template <
typename T>
835 libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
837 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
838 _coords[i] /= factor;
846 template <
typename T>
847 template <
typename T2>
853 return _coords[0]*p._coords[0];
857 return (_coords[0]*p._coords[0] +
858 _coords[1]*p._coords[1]);
862 return (_coords[0]*p(0) +
868 template <
typename T>
869 template <
typename T2>
879 template <
typename T>
880 template <
typename T2>
885 libmesh_assert_equal_to (LIBMESH_DIM, 3);
891 return TypeVector<TS>( _coords[1]*p._coords[2] - _coords[2]*p._coords[1],
892 -_coords[0]*p._coords[2] + _coords[2]*p._coords[0],
893 _coords[0]*p._coords[1] - _coords[1]*p._coords[0]);
898 #ifdef LIBMESH_ENABLE_DEPRECATED 899 template <
typename T>
903 libmesh_deprecated();
910 template <
typename T>
914 return std::sqrt(this->
norm_sq());
919 template <
typename T>
923 for (
unsigned int i=0; i<LIBMESH_DIM; i++)
929 #ifdef LIBMESH_ENABLE_DEPRECATED 930 template <
typename T>
934 libmesh_deprecated();
941 template <
typename T>
963 template <
typename T>
968 return (
std::abs(_coords[0] - rhs._coords[0])
973 return (
std::abs(_coords[0] - rhs._coords[0]) +
974 std::abs(_coords[1] - rhs._coords[1])
979 return (
std::abs(_coords[0] - rhs._coords[0]) +
980 std::abs(_coords[1] - rhs._coords[1]) +
981 std::abs(_coords[2] - rhs._coords[2])
988 template <
typename T>
993 return this->absolute_fuzzy_equals(rhs, tol *
998 return this->absolute_fuzzy_equals(rhs, tol *
1003 #if LIBMESH_DIM == 3 1004 return this->absolute_fuzzy_equals(rhs, tol *
1013 template <
typename T>
1017 #if LIBMESH_DIM == 1 1018 return (_coords[0] == rhs._coords[0]);
1021 #if LIBMESH_DIM == 2 1022 return (_coords[0] == rhs._coords[0] &&
1023 _coords[1] == rhs._coords[1]);
1026 #if LIBMESH_DIM == 3 1027 return (_coords[0] == rhs._coords[0] &&
1028 _coords[1] == rhs._coords[1] &&
1029 _coords[2] == rhs._coords[2]);
1035 template <
typename T>
1039 return (!(*
this == rhs));
1052 template <
typename T>
1058 #if LIBMESH_DIM == 3 1060 a(0)*(b(1)*c(2) - b(2)*c(1)) -
1061 a(1)*(b(0)*c(2) - b(2)*c(0)) +
1062 a(2)*(b(0)*c(1) - b(1)*c(0));
1074 template <
typename T>
1079 T z = b(0)*c(1) - b(1)*c(0);
1081 #if LIBMESH_DIM == 3 1082 T x = b(1)*c(2) - b(2)*c(1),
1083 y = b(0)*c(2) - b(2)*c(0);
1084 return x*x + y*y + z*z;
1095 template <
typename T>
1105 #endif // LIBMESH_TYPE_VECTOR_H
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator*(const Scalar) const
TypeVector< typename CompareTypes< T, T2 >::supertype > operator+(const TypeVector< T2 > &) const
CompareTypes< T, T2 >::supertype contract(const TypeVector< T2 > &) const
void add_scaled(const TypeVector< T2 > &, const T)
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
bool operator>(const TypeVector< T > &rhs) const
static const Real TOLERANCE
const TypeVector< T > & operator+=(const TypeVector< T2 > &)
void add(const TypeVector< T2 > &)
const TypeVector< T > & operator/=(const T)
T cross_norm_sq(const TypeVector< T > &b, const TypeVector< T > &c)
bool operator>=(const TypeVector< T > &rhs) const
void print(std::ostream &os=libMesh::out) const
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
TypeVector< T > unit() const
TypeVector< T > operator-() const
void write_unformatted(std::ostream &out, const bool newline=true) const
void subtract(const TypeVector< T2 > &)
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeNTensor< N, typename CompareTypes< Scalar, T >::supertype > >::type operator*(const Scalar &, const TypeNTensor< N, T > &)
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
void subtract_scaled(const TypeVector< T2 > &, const T)
const TypeVector< T > & operator-=(const TypeVector< T2 > &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void assign(const TypeVector< T2 > &)
T & slice(const unsigned int i)
const T & operator()(const unsigned int i) const
bool operator==(const TypeVector< T > &rhs) const
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator/(const Scalar) const
const TypeVector< T > & operator*=(const T)
OStreamProxy out(std::cout)
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector & >::type operator=(const Scalar &libmesh_dbg_var(p))
bool operator!=(const TypeVector< T > &rhs) const
const T & slice(const unsigned int i) const