type_vector.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_TYPE_VECTOR_H
21 #define LIBMESH_TYPE_VECTOR_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/compare_types.h"
26 #include "libmesh/tensor_tools.h"
27 
28 // C++ includes
29 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
30 #include <cmath>
31 #include <complex>
32 
33 namespace libMesh
34 {
35 
36 // Forward declarations
37 template <typename T> class TypeTensor;
38 template <typename T> class VectorValue;
39 template <typename T> class TensorValue;
40 
51 template <typename T>
52 class TypeVector
53 {
54  template <typename T2>
55  friend class TypeVector;
56 
57  friend class TypeTensor<T>;
58 
59 protected:
60 
65  TypeVector ();
66 
71  TypeVector (const T x,
72  const T y=0,
73  const T z=0);
74 
79  template <typename Scalar1, typename Scalar2, typename Scalar3>
80  TypeVector (typename
82  const Scalar1>::type x,
83  typename
85  const Scalar2>::type y=0,
86  typename
88  const Scalar3>::type z=0);
89 
96  template <typename Scalar>
97  TypeVector (const Scalar x,
98  typename
100  const Scalar>::type * sfinae = libmesh_nullptr);
101 
102 public:
103 
107  typedef T value_type;
108 
112  template <typename T2>
113  TypeVector (const TypeVector<T2> & p);
114 
118  ~TypeVector ();
119 
123  template <typename T2>
124  void assign (const TypeVector<T2> &);
125 
129  template <typename Scalar>
130  typename boostcopy::enable_if_c<
132  TypeVector &>::type
133  operator = (const Scalar & libmesh_dbg_var(p))
134  { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
135 
139  const T & operator () (const unsigned int i) const;
140  const T & slice (const unsigned int i) const { return (*this)(i); }
141 
145  T & operator () (const unsigned int i);
146  T & slice (const unsigned int i) { return (*this)(i); }
147 
153  template <typename T2>
155  operator + (const TypeVector<T2> &) const;
156 
162  template <typename T2>
163  const TypeVector<T> & operator += (const TypeVector<T2> &);
164 
168  template <typename T2>
169  void add (const TypeVector<T2> &);
170 
174  template <typename T2>
175  void add_scaled (const TypeVector<T2> &, const T);
176 
182  template <typename T2>
184  operator - (const TypeVector<T2> &) const;
185 
191  template <typename T2>
192  const TypeVector<T> & operator -= (const TypeVector<T2> &);
193 
197  template <typename T2>
198  void subtract (const TypeVector<T2> &);
199 
204  template <typename T2>
205  void subtract_scaled (const TypeVector<T2> &, const T);
206 
210  TypeVector<T> operator - () const;
211 
217  template <typename Scalar>
218  typename boostcopy::enable_if_c<
219  ScalarTraits<Scalar>::value,
221  operator * (const Scalar) const;
222 
228  const TypeVector<T> & operator *= (const T);
229 
235  template <typename Scalar>
236  typename boostcopy::enable_if_c<
237  ScalarTraits<Scalar>::value,
239  operator / (const Scalar) const;
240 
246  const TypeVector<T> & operator /= (const T);
247 
254  template <typename T2>
256  operator * (const TypeVector<T2> &) const;
257 
261  template <typename T2>
263  contract (const TypeVector<T2> &) const;
264 
268  template <typename T2>
270  cross(const TypeVector<T2> & v) const;
271 
275  TypeVector<T> unit() const;
276 
283  Real size() const;
284 
289  Real norm() const;
290 
297  Real size_sq() const;
298 
303  Real norm_sq() const;
304 
308  void zero();
309 
314  bool relative_fuzzy_equals(const TypeVector<T> & rhs,
315  Real tol = TOLERANCE) const;
316 
321  bool absolute_fuzzy_equals(const TypeVector<T> & rhs,
322  Real tol = TOLERANCE) const;
323 
331  bool operator == (const TypeVector<T> & rhs) const;
332 
336  bool operator != (const TypeVector<T> & rhs) const;
337 
344  bool operator < (const TypeVector<T> & rhs) const;
345 
352  bool operator <= (const TypeVector<T> & rhs) const;
353 
360  bool operator > (const TypeVector<T> & rhs) const;
361 
368  bool operator >= (const TypeVector<T> & rhs) const;
369 
373  void print(std::ostream & os = libMesh::out) const;
374 
383  friend std::ostream & operator << (std::ostream & os, const TypeVector<T> & t)
384  {
385  t.print(os);
386  return os;
387  }
388 
395  void write_unformatted (std::ostream & out, const bool newline = true) const;
396 
397 protected:
398 
402  T _coords[LIBMESH_DIM];
403 };
404 
405 
406 
407 
408 
409 //------------------------------------------------------
410 // Inline functions
411 
412 template <typename T>
413 inline
415 {
416  _coords[0] = 0;
417 
418 #if LIBMESH_DIM > 1
419  _coords[1] = 0;
420 #endif
421 
422 #if LIBMESH_DIM > 2
423  _coords[2] = 0;
424 #endif
425 }
426 
427 
428 
429 template <typename T>
430 inline
432  const T y,
433  const T z)
434 {
435  _coords[0] = x;
436 
437 #if LIBMESH_DIM > 1
438  _coords[1] = y;
439 #else
440  libmesh_assert_equal_to (y, 0);
441 #endif
442 
443 #if LIBMESH_DIM > 2
444  _coords[2] = z;
445 #else
446  libmesh_assert_equal_to (z, 0);
447 #endif
448 }
449 
450 
451 template <typename T>
452 template <typename Scalar1, typename Scalar2, typename Scalar3>
453 inline
456  const Scalar1>::type x,
457  typename
459  const Scalar2>::type y,
460  typename
462  const Scalar3>::type z)
463 {
464  _coords[0] = x;
465 
466 #if LIBMESH_DIM > 1
467  _coords[1] = y;
468 #else
469  libmesh_assert_equal_to (y, 0);
470 #endif
471 
472 #if LIBMESH_DIM > 2
473  _coords[2] = z;
474 #else
475  libmesh_assert_equal_to (z, 0);
476 #endif
477 }
478 
479 
480 
481 template <typename T>
482 template <typename Scalar>
483 inline
485  typename
486  boostcopy::enable_if_c<ScalarTraits<Scalar>::value,
487  const Scalar>::type * /*sfinae*/)
488 {
489  _coords[0] = x;
490 
491 #if LIBMESH_DIM > 1
492  _coords[1] = 0;
493 #endif
494 
495 #if LIBMESH_DIM > 2
496  _coords[2] = 0;
497 #endif
498 }
499 
500 
501 
502 template <typename T>
503 template <typename T2>
504 inline
506 {
507  // copy the nodes from vector p to me
508  for (unsigned int i=0; i<LIBMESH_DIM; i++)
509  _coords[i] = p._coords[i];
510 }
511 
512 
513 
514 template <typename T>
515 inline
517 {
518 }
519 
520 
521 
522 template <typename T>
523 template <typename T2>
524 inline
526 {
527  for (unsigned int i=0; i<LIBMESH_DIM; i++)
528  _coords[i] = p._coords[i];
529 }
530 
531 
532 
533 template <typename T>
534 inline
535 const T & TypeVector<T>::operator () (const unsigned int i) const
536 {
537  libmesh_assert_less (i, LIBMESH_DIM);
538 
539  return _coords[i];
540 }
541 
542 
543 
544 template <typename T>
545 inline
546 T & TypeVector<T>::operator () (const unsigned int i)
547 {
548  libmesh_assert_less (i, LIBMESH_DIM);
549 
550  return _coords[i];
551 }
552 
553 
554 
555 template <typename T>
556 template <typename T2>
557 inline
560 {
561  typedef typename CompareTypes<T, T2>::supertype TS;
562 #if LIBMESH_DIM == 1
563  return TypeVector<TS> (_coords[0] + p._coords[0]);
564 #endif
565 
566 #if LIBMESH_DIM == 2
567  return TypeVector<TS> (_coords[0] + p._coords[0],
568  _coords[1] + p._coords[1]);
569 #endif
570 
571 #if LIBMESH_DIM == 3
572  return TypeVector<TS> (_coords[0] + p._coords[0],
573  _coords[1] + p._coords[1],
574  _coords[2] + p._coords[2]);
575 #endif
576 
577 }
578 
579 
580 
581 template <typename T>
582 template <typename T2>
583 inline
585 {
586  this->add (p);
587 
588  return *this;
589 }
590 
591 
592 
593 template <typename T>
594 template <typename T2>
595 inline
597 {
598 #if LIBMESH_DIM == 1
599  _coords[0] += p._coords[0];
600 #endif
601 
602 #if LIBMESH_DIM == 2
603  _coords[0] += p._coords[0];
604  _coords[1] += p._coords[1];
605 #endif
606 
607 #if LIBMESH_DIM == 3
608  _coords[0] += p._coords[0];
609  _coords[1] += p._coords[1];
610  _coords[2] += p._coords[2];
611 #endif
612 
613 }
614 
615 
616 
617 template <typename T>
618 template <typename T2>
619 inline
620 void TypeVector<T>::add_scaled (const TypeVector<T2> & p, const T factor)
621 {
622 #if LIBMESH_DIM == 1
623  _coords[0] += factor*p(0);
624 #endif
625 
626 #if LIBMESH_DIM == 2
627  _coords[0] += factor*p(0);
628  _coords[1] += factor*p(1);
629 #endif
630 
631 #if LIBMESH_DIM == 3
632  _coords[0] += factor*p(0);
633  _coords[1] += factor*p(1);
634  _coords[2] += factor*p(2);
635 #endif
636 
637 }
638 
639 
640 
641 template <typename T>
642 template <typename T2>
643 inline
646 {
647  typedef typename CompareTypes<T, T2>::supertype TS;
648 
649 #if LIBMESH_DIM == 1
650  return TypeVector<TS>(_coords[0] - p._coords[0]);
651 #endif
652 
653 #if LIBMESH_DIM == 2
654  return TypeVector<TS>(_coords[0] - p._coords[0],
655  _coords[1] - p._coords[1]);
656 #endif
657 
658 #if LIBMESH_DIM == 3
659  return TypeVector<TS>(_coords[0] - p._coords[0],
660  _coords[1] - p._coords[1],
661  _coords[2] - p._coords[2]);
662 #endif
663 
664 }
665 
666 
667 
668 template <typename T>
669 template <typename T2>
670 inline
672 {
673  this->subtract (p);
674 
675  return *this;
676 }
677 
678 
679 
680 template <typename T>
681 template <typename T2>
682 inline
684 {
685  for (unsigned int i=0; i<LIBMESH_DIM; i++)
686  _coords[i] -= p._coords[i];
687 }
688 
689 
690 
691 template <typename T>
692 template <typename T2>
693 inline
694 void TypeVector<T>::subtract_scaled (const TypeVector<T2> & p, const T factor)
695 {
696  for (unsigned int i=0; i<LIBMESH_DIM; i++)
697  _coords[i] -= factor*p(i);
698 }
699 
700 
701 
702 template <typename T>
703 inline
705 {
706 
707 #if LIBMESH_DIM == 1
708  return TypeVector(-_coords[0]);
709 #endif
710 
711 #if LIBMESH_DIM == 2
712  return TypeVector(-_coords[0],
713  -_coords[1]);
714 #endif
715 
716 #if LIBMESH_DIM == 3
717  return TypeVector(-_coords[0],
718  -_coords[1],
719  -_coords[2]);
720 #endif
721 
722 }
723 
724 
725 
726 template <typename T>
727 template <typename Scalar>
728 inline
729 typename boostcopy::enable_if_c<
730  ScalarTraits<Scalar>::value,
732 TypeVector<T>::operator * (const Scalar factor) const
733 {
734  typedef typename CompareTypes<T, Scalar>::supertype SuperType;
735 
736 #if LIBMESH_DIM == 1
737  return TypeVector<SuperType>(_coords[0]*factor);
738 #endif
739 
740 #if LIBMESH_DIM == 2
741  return TypeVector<SuperType>(_coords[0]*factor,
742  _coords[1]*factor);
743 #endif
744 
745 #if LIBMESH_DIM == 3
746  return TypeVector<SuperType>(_coords[0]*factor,
747  _coords[1]*factor,
748  _coords[2]*factor);
749 #endif
750 }
751 
752 
753 
754 template <typename T, typename Scalar>
755 inline
756 typename boostcopy::enable_if_c<
757  ScalarTraits<Scalar>::value,
759 operator * (const Scalar factor,
760  const TypeVector<T> & v)
761 {
762  return v * factor;
763 }
764 
765 
766 
767 template <typename T>
768 inline
770 {
771 #if LIBMESH_DIM == 1
772  _coords[0] *= factor;
773 #endif
774 
775 #if LIBMESH_DIM == 2
776  _coords[0] *= factor;
777  _coords[1] *= factor;
778 #endif
779 
780 #if LIBMESH_DIM == 3
781  _coords[0] *= factor;
782  _coords[1] *= factor;
783  _coords[2] *= factor;
784 #endif
785 
786  return *this;
787 }
788 
789 
790 
791 template <typename T>
792 template <typename Scalar>
793 inline
794 typename boostcopy::enable_if_c<
795  ScalarTraits<Scalar>::value,
797 TypeVector<T>::operator / (const Scalar factor) const
798 {
799  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
800 
801  typedef typename CompareTypes<T, Scalar>::supertype TS;
802 
803 #if LIBMESH_DIM == 1
804  return TypeVector<TS>(_coords[0]/factor);
805 #endif
806 
807 #if LIBMESH_DIM == 2
808  return TypeVector<TS>(_coords[0]/factor,
809  _coords[1]/factor);
810 #endif
811 
812 #if LIBMESH_DIM == 3
813  return TypeVector<TS>(_coords[0]/factor,
814  _coords[1]/factor,
815  _coords[2]/factor);
816 #endif
817 
818 }
819 
820 
821 
822 
823 template <typename T>
824 inline
825 const TypeVector<T> &
827 {
828  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
829 
830  for (unsigned int i=0; i<LIBMESH_DIM; i++)
831  _coords[i] /= factor;
832 
833  return *this;
834 }
835 
836 
837 
838 
839 template <typename T>
840 template <typename T2>
841 inline
844 {
845 #if LIBMESH_DIM == 1
846  return _coords[0]*p._coords[0];
847 #endif
848 
849 #if LIBMESH_DIM == 2
850  return (_coords[0]*p._coords[0] +
851  _coords[1]*p._coords[1]);
852 #endif
853 
854 #if LIBMESH_DIM == 3
855  return (_coords[0]*p(0) +
856  _coords[1]*p(1) +
857  _coords[2]*p(2));
858 #endif
859 }
860 
861 template <typename T>
862 template <typename T2>
863 inline
866 {
867  return (*this)*(p);
868 }
869 
870 
871 
872 template <typename T>
873 template <typename T2>
876 {
877  typedef typename CompareTypes<T, T2>::supertype TS;
878  libmesh_assert_equal_to (LIBMESH_DIM, 3);
879 
880  // | i j k |
881  // |(*this)(0) (*this)(1) (*this)(2)|
882  // | p(0) p(1) p(2) |
883 
884  return TypeVector<TS>( _coords[1]*p._coords[2] - _coords[2]*p._coords[1],
885  -_coords[0]*p._coords[2] + _coords[2]*p._coords[0],
886  _coords[0]*p._coords[1] - _coords[1]*p._coords[0]);
887 }
888 
889 
890 
891 template <typename T>
892 inline
894 {
895  libmesh_deprecated();
896  return this->norm();
897 }
898 
899 
900 
901 template <typename T>
902 inline
904 {
905  return std::sqrt(this->norm_sq());
906 }
907 
908 
909 
910 template <typename T>
911 inline
913 {
914  for (unsigned int i=0; i<LIBMESH_DIM; i++)
915  _coords[i] = 0.;
916 }
917 
918 
919 
920 template <typename T>
921 inline
923 {
924  libmesh_deprecated();
925  return this->norm_sq();
926 }
927 
928 
929 
930 template <typename T>
931 inline
933 {
934 #if LIBMESH_DIM == 1
935  return (TensorTools::norm_sq(_coords[0]));
936 #endif
937 
938 #if LIBMESH_DIM == 2
939  return (TensorTools::norm_sq(_coords[0]) +
941 #endif
942 
943 #if LIBMESH_DIM == 3
944  return (TensorTools::norm_sq(_coords[0]) +
947 #endif
948 }
949 
950 
951 
952 template <typename T>
953 inline
955 {
956 #if LIBMESH_DIM == 1
957  return (std::abs(_coords[0] - rhs._coords[0])
958  <= tol);
959 #endif
960 
961 #if LIBMESH_DIM == 2
962  return (std::abs(_coords[0] - rhs._coords[0]) +
963  std::abs(_coords[1] - rhs._coords[1])
964  <= tol);
965 #endif
966 
967 #if LIBMESH_DIM == 3
968  return (std::abs(_coords[0] - rhs._coords[0]) +
969  std::abs(_coords[1] - rhs._coords[1]) +
970  std::abs(_coords[2] - rhs._coords[2])
971  <= tol);
972 #endif
973 }
974 
975 
976 
977 template <typename T>
978 inline
980 {
981 #if LIBMESH_DIM == 1
982  return this->absolute_fuzzy_equals(rhs, tol *
983  (std::abs(_coords[0]) + std::abs(rhs._coords[0])));
984 #endif
985 
986 #if LIBMESH_DIM == 2
987  return this->absolute_fuzzy_equals(rhs, tol *
988  (std::abs(_coords[0]) + std::abs(rhs._coords[0]) +
989  std::abs(_coords[1]) + std::abs(rhs._coords[1])));
990 #endif
991 
992 #if LIBMESH_DIM == 3
993  return this->absolute_fuzzy_equals(rhs, tol *
994  (std::abs(_coords[0]) + std::abs(rhs._coords[0]) +
995  std::abs(_coords[1]) + std::abs(rhs._coords[1]) +
996  std::abs(_coords[2]) + std::abs(rhs._coords[2])));
997 #endif
998 }
999 
1000 
1001 
1002 template <typename T>
1003 inline
1005 {
1006 #if LIBMESH_DIM == 1
1007  return (_coords[0] == rhs._coords[0]);
1008 #endif
1009 
1010 #if LIBMESH_DIM == 2
1011  return (_coords[0] == rhs._coords[0] &&
1012  _coords[1] == rhs._coords[1]);
1013 #endif
1014 
1015 #if LIBMESH_DIM == 3
1016  return (_coords[0] == rhs._coords[0] &&
1017  _coords[1] == rhs._coords[1] &&
1018  _coords[2] == rhs._coords[2]);
1019 #endif
1020 }
1021 
1022 
1023 
1024 template <typename T>
1025 inline
1027 {
1028  return (!(*this == rhs));
1029 }
1030 
1031 
1032 //------------------------------------------------------
1033 // Non-member functions on TypeVectors
1034 
1035 // Compute a * (b.cross(c)) without creating a temporary
1036 // for the cross product. Equivalent to the determinant
1037 // of the 3x3 tensor:
1038 // [a0, a1, a2]
1039 // [b0, b1, b2]
1040 // [c0, c1, c2]
1041 template <typename T>
1042 inline
1044  const TypeVector<T> & b,
1045  const TypeVector<T> & c)
1046 {
1047  // We only support cross products when LIBMESH_DIM==3, same goes for this.
1048  libmesh_assert_equal_to (LIBMESH_DIM, 3);
1049 
1050  return
1051  a(0)*(b(1)*c(2) - b(2)*c(1)) -
1052  a(1)*(b(0)*c(2) - b(2)*c(0)) +
1053  a(2)*(b(0)*c(1) - b(1)*c(0));
1054 }
1055 
1056 
1057 
1062 template <typename T>
1063 inline
1065  const TypeVector<T> & c)
1066 {
1067  // We only support cross products when LIBMESH_DIM==3, same goes for this.
1068  libmesh_assert_equal_to (LIBMESH_DIM, 3);
1069 
1070  T x = b(1)*c(2) - b(2)*c(1),
1071  y = b(0)*c(2) - b(2)*c(0),
1072  z = b(0)*c(1) - b(1)*c(0);
1073 
1074  return x*x + y*y + z*z;
1075 }
1076 
1077 
1078 
1082 template <typename T>
1083 inline
1085  const TypeVector<T> & c)
1086 {
1087  // We only support cross products when LIBMESH_DIM==3, same goes for this.
1088  libmesh_assert_equal_to (LIBMESH_DIM, 3);
1089 
1090  return std::sqrt(cross_norm_sq(b,c));
1091 }
1092 
1093 } // namespace libMesh
1094 
1095 #endif // LIBMESH_TYPE_VECTOR_H
bool operator>=(const TypeVector< T > &rhs) const
Definition: type_vector.C:152
T _coords[LIBMESH_DIM]
Definition: type_vector.h:402
double abs(double a)
bool operator==(const TypeVector< T > &rhs) const
Definition: type_vector.h:1004
const T & slice(const unsigned int i) const
Definition: type_vector.h:140
bool operator>(const TypeVector< T > &rhs) const
Definition: type_vector.C:138
CompareTypes< T, T2 >::supertype contract(const TypeVector< T2 > &) const
Definition: type_vector.h:865
void add_scaled(const TypeVector< T2 > &, const T)
Definition: type_vector.h:620
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1084
void write_unformatted(std::ostream &out, const bool newline=true) const
Definition: type_vector.C:92
Real norm() const
Definition: type_vector.h:903
const class libmesh_nullptr_t libmesh_nullptr
static const Real TOLERANCE
TypeVector< T > operator-() const
Definition: type_vector.h:704
TypeVector< typename CompareTypes< T, T2 >::supertype > operator+(const TypeVector< T2 > &) const
Definition: type_vector.h:559
const TypeVector< T > & operator+=(const TypeVector< T2 > &)
Definition: type_vector.h:584
void add(const TypeVector< T2 > &)
Definition: type_vector.h:596
const TypeVector< T > & operator/=(const T)
Definition: type_vector.h:826
T cross_norm_sq(const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1064
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1043
PetscErrorCode Vec x
void subtract(const TypeVector< T2 > &)
Definition: type_vector.h:683
void print(std::ostream &os=libMesh::out) const
Definition: type_vector.C:64
bool operator!=(const TypeVector< T > &rhs) const
Definition: type_vector.h:1026
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:875
PetscErrorCode Vec Mat libmesh_dbg_var(j)
Real norm_sq() const
Definition: type_vector.h:932
void subtract_scaled(const TypeVector< T2 > &, const T)
Definition: type_vector.h:694
TypeVector< T > unit() const
Definition: type_vector.C:37
const TypeVector< T > & operator-=(const TypeVector< T2 > &)
Definition: type_vector.h:671
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:954
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator/(const Scalar) const
Definition: type_vector.h:797
void assign(const TypeVector< T2 > &)
Definition: type_vector.h:525
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator*(const Scalar) const
Definition: type_vector.h:732
const T & operator()(const unsigned int i) const
Definition: type_vector.h:535
T & slice(const unsigned int i)
Definition: type_vector.h:146
Real size_sq() const
Definition: type_vector.h:922
Real size() const
Definition: type_vector.h:893
const TypeVector< T > & operator*=(const T)
Definition: type_vector.h:769
OStreamProxy out(std::cout)
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector & >::type operator=(const Scalar &libmesh_dbg_var(p))
Definition: type_vector.h:133
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:979