libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
type_vector.h
浏览该文件的文档.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2023 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 #ifdef LIBMESH_HAVE_EIGEN
34 #include "libmesh/ignore_warnings.h"
35 #include "Eigen/Core"
36 #include "libmesh/restore_warnings.h"
37 #endif
38 
39 #ifdef LIBMESH_HAVE_METAPHYSICL
40 #include "metaphysicl/dualnumber_decl.h"
41 
42 namespace std
43 {
44 #ifdef LIBMESH_HAVE_EIGEN
45 
51 template <typename T, typename D, bool asd>
52 using ADRealEigenVector = Eigen::Matrix<MetaPhysicL::DualNumber<T, D, asd>, Eigen::Dynamic, 1>;
53 
63 template <typename T, typename D, bool asd>
64 ADRealEigenVector<T, D, asd> norm(const ADRealEigenVector<T, D, asd> &) {throw std::runtime_error("unimplemented");}
65 
75 template <typename T, typename D, bool asd>
76 ADRealEigenVector<T, D, asd> norm(ADRealEigenVector<T, D, asd> &&) {throw std::runtime_error("unimplemented");}
77 
87 template <typename T, typename D, bool asd>
88 ADRealEigenVector<T, D, asd> sqrt(const ADRealEigenVector<T, D, asd> &) {throw std::runtime_error("unimplemented");}
89 
99 template <typename T, typename D, bool asd>
100 ADRealEigenVector<T, D, asd> sqrt(ADRealEigenVector<T, D, asd> &&) {throw std::runtime_error("unimplemented");}
101 
111 template <typename T, typename D, bool asd>
112 ADRealEigenVector<T, D, asd> abs(const ADRealEigenVector<T, D, asd> &) {throw std::runtime_error("unimplemented");}
113 
123 template <typename T, typename D, bool asd>
124 ADRealEigenVector<T, D, asd> abs(ADRealEigenVector<T, D, asd> &&) {throw std::runtime_error("unimplemented");}
125 #endif
126 }
127 #endif
128 
129 namespace libMesh
130 {
131 // Forward declarations
132 template <typename T> class TypeTensor;
133 template <typename T> class VectorValue;
134 template <typename T> class TensorValue;
135 }
136 
137 namespace std
138 {
139 template <typename T>
140 auto norm(const libMesh::TypeVector<T> & vector) -> decltype(std::norm(T()));
141 }
142 
143 namespace libMesh
144 {
153 template <typename T>
154 class TypeVector
155 {
156  template <typename T2>
157  friend class TypeVector;
158 
159  friend class TypeTensor<T>;
160 
161 public:
165  TypeVector ();
166 
167 protected:
171  TypeVector (const T & x,
172  const T & y=0,
173  const T & z=0);
174 
178  template <typename Scalar1, typename Scalar2, typename Scalar3>
179  TypeVector (typename
180  boostcopy::enable_if_c<ScalarTraits<Scalar1>::value,
181  const Scalar1>::type & x,
182  typename
183  boostcopy::enable_if_c<ScalarTraits<Scalar2>::value,
184  const Scalar2>::type & y=0,
185  typename
186  boostcopy::enable_if_c<ScalarTraits<Scalar3>::value,
187  const Scalar3>::type & z=0);
188 
195  template <typename Scalar>
196  TypeVector (const Scalar & x,
197  typename
198  boostcopy::enable_if_c<ScalarTraits<Scalar>::value,
199  const Scalar>::type * sfinae = nullptr);
200 
201 public:
202 
206  typedef T value_type;
207 
211  typedef unsigned int index_type;
212 
219  template <typename T2>
220  TypeVector (const TypeVector<T2> & p);
221 
225  TypeVector (const TypeVector & p) = default;
226 
230  ~TypeVector () = default;
231 
238  template <typename T2>
239  void assign (const TypeVector<T2> &);
240 
248  template <typename Scalar>
249  typename boostcopy::enable_if_c<
250  ScalarTraits<Scalar>::value,
251  TypeVector &>::type
252  operator = (const Scalar & libmesh_dbg_var(p))
253  { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
254 
261  const T & operator () (const unsigned int i) const;
262  const T & slice (const unsigned int i) const { return (*this)(i); }
263 
270  T & operator () (const unsigned int i);
271  T & slice (const unsigned int i) { return (*this)(i); }
272 
280  template <typename T2>
282  operator + (const TypeVector<T2> &) const;
283 
291  template <typename T2>
292  const TypeVector<T> & operator += (const TypeVector<T2> &);
293 
301  template <typename T2>
303  operator-(const TypeVector<T2> &) const;
304 
312  template <typename T2>
313  const TypeVector<T> &operator-=(const TypeVector<T2> &);
314 
321  template <typename T2>
322  void subtract(const TypeVector<T2> &);
323 
331  template <typename T2>
332  void subtract_scaled(const TypeVector<T2> &, const T &);
333 
339  TypeVector<T> operator-() const;
340 
348  template <typename Scalar>
349  typename boostcopy::enable_if_c<
350  ScalarTraits<Scalar>::value,
352  operator*(const Scalar &) const;
353 
361  const TypeVector<T> &operator*=(const T &);
362 
370  template <typename Scalar>
371  typename boostcopy::enable_if_c<
372  ScalarTraits<Scalar>::value,
374  operator/(const Scalar &) const;
375 
383  const TypeVector<T> &operator/=(const T &);
384 
395 template <typename T2>
396 typename CompareTypes<T, T2>::supertype
397 operator*(const TypeVector<T2> &) const;
398 
406 template <typename T2>
407 typename CompareTypes<T, T2>::supertype
408 contract(const TypeVector<T2> &) const;
409 
417 template <typename T2>
419 cross(const TypeVector<T2> &v) const;
420 
426 TypeVector<T> unit() const;
427 
433 auto norm() const -> decltype(std::norm(T()));
434 
440 auto norm_sq() const -> decltype(std::norm(T()));
441 
447 bool is_zero() const;
448 
452 void zero();
453 
461 bool relative_fuzzy_equals(const TypeVector<T> &rhs, Real tol = TOLERANCE) const;
462 
470 bool absolute_fuzzy_equals(const TypeVector<T> &rhs, Real tol = TOLERANCE) const;
471 
480 bool operator==(const TypeVector<T> &rhs) const;
481 
488 bool operator!=(const TypeVector<T> &rhs) const;
489 
498 bool operator<(const TypeVector<T> &rhs) const;
499 
508 bool operator<=(const TypeVector<T> &rhs) const;
509 
518 bool operator>(const TypeVector<T> &rhs) const;
519 
528 bool operator>=(const TypeVector<T> &rhs) const;
529 
535 void print(std::ostream &os = libMesh::out) const;
536 
545 friend std::ostream &operator<<(std::ostream &os, const TypeVector<T> &t)
546 {
547  t.print(os);
548  return os;
549 }
550 
558 void write_unformatted(std::ostream &out_stream, const bool newline = true) const;
559 
560 protected:
561 
565 T _coords[LIBMESH_DIM];
566 };
567 
568 
569 
570 
571 
572 //------------------------------------------------------
573 // Inline functions
574 
575 template <typename T>
576 inline
578 {
579  _coords[0] = {};
580 
581 #if LIBMESH_DIM > 1
582  _coords[1] = {};
583 #endif
584 
585 #if LIBMESH_DIM > 2
586  _coords[2] = {};
587 #endif
588 }
589 
590 
591 
592 template <typename T>
593 inline
595  const T & y,
596  const T & z)
597 {
598  _coords[0] = x;
599 
600 #if LIBMESH_DIM > 1
601  _coords[1] = y;
602 #else
603  libmesh_ignore(y);
604  libmesh_assert_equal_to (y, 0);
605 #endif
606 
607 #if LIBMESH_DIM > 2
608  _coords[2] = z;
609 #else
610  libmesh_ignore(z);
611  libmesh_assert_equal_to (z, 0);
612 #endif
613 }
614 
615 
616 template <typename T>
617 template <typename Scalar1, typename Scalar2, typename Scalar3>
618 inline
620  boostcopy::enable_if_c<ScalarTraits<Scalar1>::value,
621  const Scalar1>::type & x,
622  typename
623  boostcopy::enable_if_c<ScalarTraits<Scalar2>::value,
624  const Scalar2>::type & y,
625  typename
626  boostcopy::enable_if_c<ScalarTraits<Scalar3>::value,
627  const Scalar3>::type & z)
628 {
629  _coords[0] = x;
630 
631 #if LIBMESH_DIM > 1
632  _coords[1] = y;
633 #else
634  libmesh_assert_equal_to (y, 0);
635 #endif
636 
637 #if LIBMESH_DIM > 2
638  _coords[2] = z;
639 #else
640  libmesh_assert_equal_to (z, 0);
641 #endif
642 }
643 
644 
645 
646 template <typename T>
647 template <typename Scalar>
648 inline
649 TypeVector<T>::TypeVector (const Scalar & x,
650  typename
651  boostcopy::enable_if_c<ScalarTraits<Scalar>::value,
652  const Scalar>::type * /*sfinae*/)
653 {
654  _coords[0] = x;
655 
656 #if LIBMESH_DIM > 1
657  _coords[1] = 0;
658 #endif
659 
660 #if LIBMESH_DIM > 2
661  _coords[2] = 0;
662 #endif
663 }
664 
665 
666 
667 template <typename T>
668 template <typename T2>
669 inline
671 {
672  // copy the nodes from vector p to me
673  for (unsigned int i=0; i<LIBMESH_DIM; i++)
674  _coords[i] = p._coords[i];
675 }
676 
677 
678 
679 template <typename T>
680 template <typename T2>
681 inline
683 {
684  for (unsigned int i=0; i<LIBMESH_DIM; i++)
685  _coords[i] = p._coords[i];
686 }
687 
688 
689 
690 template <typename T>
691 inline
692 const T & TypeVector<T>::operator () (const unsigned int i) const
693 {
694  libmesh_assert_less (i, LIBMESH_DIM);
695 
696  return _coords[i];
697 }
698 
699 
700 
701 template <typename T>
702 inline
703 T & TypeVector<T>::operator () (const unsigned int i)
704 {
705  libmesh_assert_less (i, LIBMESH_DIM);
706 
707  return _coords[i];
708 }
709 
710 
711 
712 template <typename T>
713 template <typename T2>
714 inline
717 {
718  typedef typename CompareTypes<T, T2>::supertype TS;
719 #if LIBMESH_DIM == 1
720  return TypeVector<TS> (_coords[0] + p._coords[0]);
721 #endif
722 
723 #if LIBMESH_DIM == 2
724  return TypeVector<TS> (_coords[0] + p._coords[0],
725  _coords[1] + p._coords[1]);
726 #endif
727 
728 #if LIBMESH_DIM == 3
729  return TypeVector<TS> (_coords[0] + p._coords[0],
730  _coords[1] + p._coords[1],
731  _coords[2] + p._coords[2]);
732 #endif
733 
734 }
735 
736 
737 
738 template <typename T>
739 template <typename T2>
740 inline
742 {
743  this->add (p);
744 
745  return *this;
746 }
747 
748 
749 
750 template <typename T>
751 template <typename T2>
752 inline
753 void TypeVector<T>::add (const TypeVector<T2> & p)
754 {
755 #if LIBMESH_DIM == 1
756  _coords[0] += p._coords[0];
757 #endif
758 
759 #if LIBMESH_DIM == 2
760  _coords[0] += p._coords[0];
761  _coords[1] += p._coords[1];
762 #endif
763 
764 #if LIBMESH_DIM == 3
765  _coords[0] += p._coords[0];
766  _coords[1] += p._coords[1];
767  _coords[2] += p._coords[2];
768 #endif
769 
770 }
771 
772 
773 
774 template <typename T>
775 template <typename T2>
776 inline
777 void TypeVector<T>::add_scaled (const TypeVector<T2> & p, const T & factor)
778 {
779 #if LIBMESH_DIM == 1
780  _coords[0] += factor*p(0);
781 #endif
782 
783 #if LIBMESH_DIM == 2
784  _coords[0] += factor*p(0);
785  _coords[1] += factor*p(1);
786 #endif
787 
788 #if LIBMESH_DIM == 3
789  _coords[0] += factor*p(0);
790  _coords[1] += factor*p(1);
791  _coords[2] += factor*p(2);
792 #endif
793 
794 }
795 
796 
797 
798 template <typename T>
799 template <typename T2>
800 inline
801 TypeVector<typename CompareTypes<T, T2>::supertype>
803 {
804  typedef typename CompareTypes<T, T2>::supertype TS;
805 
806 #if LIBMESH_DIM == 1
807  return TypeVector<TS>(_coords[0] - p._coords[0]);
808 #endif
809 
810 #if LIBMESH_DIM == 2
811  return TypeVector<TS>(_coords[0] - p._coords[0],
812  _coords[1] - p._coords[1]);
813 #endif
814 
815 #if LIBMESH_DIM == 3
816  return TypeVector<TS>(_coords[0] - p._coords[0],
817  _coords[1] - p._coords[1],
818  _coords[2] - p._coords[2]);
819 #endif
820 
821 }
822 
823 
824 
825 template <typename T>
826 template <typename T2>
827 inline
829 {
830  this->subtract (p);
831 
832  return *this;
833 }
834 
835 
836 
837 template <typename T>
838 template <typename T2>
839 inline
841 {
842  for (unsigned int i=0; i<LIBMESH_DIM; i++)
843  _coords[i] -= p._coords[i];
844 }
845 
846 
847 
848 template <typename T>
849 template <typename T2>
850 inline
851 void TypeVector<T>::subtract_scaled (const TypeVector<T2> & p, const T & factor)
852 {
853  for (unsigned int i=0; i<LIBMESH_DIM; i++)
854  _coords[i] -= factor*p(i);
855 }
856 
857 
858 
859 template <typename T>
860 inline
862 {
863 
864 #if LIBMESH_DIM == 1
865  return TypeVector(-_coords[0]);
866 #endif
867 
868 #if LIBMESH_DIM == 2
869  return TypeVector(-_coords[0],
870  -_coords[1]);
871 #endif
872 
873 #if LIBMESH_DIM == 3
874  return TypeVector(-_coords[0],
875  -_coords[1],
876  -_coords[2]);
877 #endif
878 
879 }
880 
881 
882 
883 template <typename T>
884 template <typename Scalar>
885 inline
886 typename boostcopy::enable_if_c<
887  ScalarTraits<Scalar>::value,
889 TypeVector<T>::operator * (const Scalar & factor) const
890 {
891  typedef typename CompareTypes<T, Scalar>::supertype SuperType;
892 
893 #if LIBMESH_DIM == 1
894  return TypeVector<SuperType>(_coords[0]*factor);
895 #endif
896 
897 #if LIBMESH_DIM == 2
898  return TypeVector<SuperType>(_coords[0]*factor,
899  _coords[1]*factor);
900 #endif
901 
902 #if LIBMESH_DIM == 3
903  return TypeVector<SuperType>(_coords[0]*factor,
904  _coords[1]*factor,
905  _coords[2]*factor);
906 #endif
907 }
908 
909 
910 
911 template <typename T, typename Scalar>
912 inline
913 typename boostcopy::enable_if_c<
914  ScalarTraits<Scalar>::value,
916 operator * (const Scalar & factor,
917  const TypeVector<T> & v)
918 {
919  return v * factor;
920 }
921 
922 
923 
924 template <typename T>
925 inline
927 {
928 #if LIBMESH_DIM == 1
929  _coords[0] *= factor;
930 #endif
931 
932 #if LIBMESH_DIM == 2
933  _coords[0] *= factor;
934  _coords[1] *= factor;
935 #endif
936 
937 #if LIBMESH_DIM == 3
938  _coords[0] *= factor;
939  _coords[1] *= factor;
940  _coords[2] *= factor;
941 #endif
942 
943  return *this;
944 }
945 
946 
947 
948 template <typename T>
949 template <typename Scalar>
950 inline
951 typename boostcopy::enable_if_c<
952  ScalarTraits<Scalar>::value,
954 TypeVector<T>::operator / (const Scalar & factor) const
955 {
956  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
957 
958  typedef typename CompareTypes<T, Scalar>::supertype TS;
959 
960 #if LIBMESH_DIM == 1
961  return TypeVector<TS>(_coords[0]/factor);
962 #endif
963 
964 #if LIBMESH_DIM == 2
965  return TypeVector<TS>(_coords[0]/factor,
966  _coords[1]/factor);
967 #endif
968 
969 #if LIBMESH_DIM == 3
970  return TypeVector<TS>(_coords[0]/factor,
971  _coords[1]/factor,
972  _coords[2]/factor);
973 #endif
974 
975 }
976 
977 
978 
979 
980 template <typename T>
981 inline
982 const TypeVector<T> &
984 {
985  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
986 
987  for (unsigned int i=0; i<LIBMESH_DIM; i++)
988  _coords[i] /= factor;
989 
990  return *this;
991 }
992 
993 
994 
995 
996 template <typename T>
997 template <typename T2>
998 inline
999 typename CompareTypes<T, T2>::supertype
1001 {
1002 #if LIBMESH_DIM == 1
1003  return _coords[0]*p._coords[0];
1004 #endif
1005 
1006 #if LIBMESH_DIM == 2
1007  return (_coords[0]*p._coords[0] +
1008  _coords[1]*p._coords[1]);
1009 #endif
1010 
1011 #if LIBMESH_DIM == 3
1012  return (_coords[0]*p(0) +
1013  _coords[1]*p(1) +
1014  _coords[2]*p(2));
1015 #endif
1016 }
1017 
1018 template <typename T>
1019 template <typename T2>
1020 inline
1021 typename CompareTypes<T, T2>::supertype
1023 {
1024  return (*this)*(p);
1025 }
1026 
1027 
1028 
1029 template <typename T>
1030 template <typename T2>
1033 {
1034  typedef typename CompareTypes<T, T2>::supertype TS;
1035  libmesh_assert_equal_to (LIBMESH_DIM, 3);
1036 
1037  // | i j k |
1038  // |(*this)(0) (*this)(1) (*this)(2)|
1039  // | p(0) p(1) p(2) |
1040 
1041 #if LIBMESH_DIM == 3
1042  return TypeVector<TS>( _coords[1]*p._coords[2] - _coords[2]*p._coords[1],
1043  -_coords[0]*p._coords[2] + _coords[2]*p._coords[0],
1044  _coords[0]*p._coords[1] - _coords[1]*p._coords[0]);
1045 #else
1046  libmesh_ignore(p);
1047  return TypeVector<TS>(0);
1048 #endif
1049 }
1050 
1051 
1052 
1053 template <typename T>
1054 inline
1055 auto TypeVector<T>::norm() const -> decltype(std::norm(T()))
1056 {
1057  return std::sqrt(this->norm_sq());
1058 }
1059 
1060 
1061 
1062 template <typename T>
1063 inline
1065 {
1066  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1067  _coords[i] = 0.;
1068 }
1069 
1070 
1071 
1072 template <typename T>
1073 inline
1074 auto TypeVector<T>::norm_sq() const -> decltype(std::norm(T()))
1075 {
1076 #if LIBMESH_DIM == 1
1077  return (TensorTools::norm_sq(_coords[0]));
1078 #endif
1079 
1080 #if LIBMESH_DIM == 2
1081  return (TensorTools::norm_sq(_coords[0]) +
1082  TensorTools::norm_sq(_coords[1]));
1083 #endif
1084 
1085 #if LIBMESH_DIM == 3
1086  return (TensorTools::norm_sq(_coords[0]) +
1087  TensorTools::norm_sq(_coords[1]) +
1088  TensorTools::norm_sq(_coords[2]));
1089 #endif
1090 }
1091 
1092 
1093 template <typename T>
1094 inline
1096 {
1097  for (const auto & val : _coords)
1098  if (val != T(0))
1099  return false;
1100  return true;
1101 }
1102 
1103 template <typename T>
1104 inline
1106 {
1107 #if LIBMESH_DIM == 1
1108  return (std::abs(_coords[0] - rhs._coords[0])
1109  <= tol);
1110 #endif
1111 
1112 #if LIBMESH_DIM == 2
1113  return (std::abs(_coords[0] - rhs._coords[0]) +
1114  std::abs(_coords[1] - rhs._coords[1])
1115  <= tol);
1116 #endif
1117 
1118 #if LIBMESH_DIM == 3
1119  return (std::abs(_coords[0] - rhs._coords[0]) +
1120  std::abs(_coords[1] - rhs._coords[1]) +
1121  std::abs(_coords[2] - rhs._coords[2])
1122  <= tol);
1123 #endif
1124 }
1125 
1126 
1127 
1128 template <typename T>
1129 inline
1131 {
1132 #if LIBMESH_DIM == 1
1133  return this->absolute_fuzzy_equals(rhs, tol *
1134  (std::abs(_coords[0]) + std::abs(rhs._coords[0])));
1135 #endif
1136 
1137 #if LIBMESH_DIM == 2
1138  return this->absolute_fuzzy_equals(rhs, tol *
1139  (std::abs(_coords[0]) + std::abs(rhs._coords[0]) +
1140  std::abs(_coords[1]) + std::abs(rhs._coords[1])));
1141 #endif
1142 
1143 #if LIBMESH_DIM == 3
1144  return this->absolute_fuzzy_equals(rhs, tol *
1145  (std::abs(_coords[0]) + std::abs(rhs._coords[0]) +
1146  std::abs(_coords[1]) + std::abs(rhs._coords[1]) +
1147  std::abs(_coords[2]) + std::abs(rhs._coords[2])));
1148 #endif
1149 }
1150 
1151 
1152 
1153 template <typename T>
1154 inline
1156 {
1157 #if LIBMESH_DIM == 1
1158  return (_coords[0] == rhs._coords[0]);
1159 #endif
1160 
1161 #if LIBMESH_DIM == 2
1162  return (_coords[0] == rhs._coords[0] &&
1163  _coords[1] == rhs._coords[1]);
1164 #endif
1165 
1166 #if LIBMESH_DIM == 3
1167  return (_coords[0] == rhs._coords[0] &&
1168  _coords[1] == rhs._coords[1] &&
1169  _coords[2] == rhs._coords[2]);
1170 #endif
1171 }
1172 
1173 
1174 
1175 template <typename T>
1176 inline
1178 {
1179  return (!(*this == rhs));
1180 }
1181 
1182 
1183 //------------------------------------------------------
1184 // Non-member functions on TypeVectors
1185 
1186 // Compute a * (b.cross(c)) without creating a temporary
1187 // for the cross product. Equivalent to the determinant
1188 // of the 3x3 tensor:
1189 // [a0, a1, a2]
1190 // [b0, b1, b2]
1191 // [c0, c1, c2]
1192 template <typename T>
1193 inline
1195  const TypeVector<T> & b,
1196  const TypeVector<T> & c)
1197 {
1198 #if LIBMESH_DIM == 3
1199  return
1200  a(0)*(b(1)*c(2) - b(2)*c(1)) -
1201  a(1)*(b(0)*c(2) - b(2)*c(0)) +
1202  a(2)*(b(0)*c(1) - b(1)*c(0));
1203 #else
1204  libmesh_ignore(a, b, c);
1205  return 0;
1206 #endif
1207 }
1208 
1209 
1210 
1215 template <typename T>
1216 inline
1218  const TypeVector<T> & c)
1219 {
1220  T z = b(0)*c(1) - b(1)*c(0);
1221 
1222 #if LIBMESH_DIM == 3
1223  T x = b(1)*c(2) - b(2)*c(1),
1224  y = b(0)*c(2) - b(2)*c(0);
1225  return x*x + y*y + z*z;
1226 #else
1227  return z*z;
1228 #endif
1229 }
1230 
1231 
1232 
1236 template <typename T>
1237 inline
1239  const TypeVector<T> & c)
1240 {
1241  return std::sqrt(cross_norm_sq(b,c));
1242 }
1243 
1244 template <typename T>
1245 inline
1247 {
1248 
1249  auto && length = norm();
1250 
1251  libmesh_assert_not_equal_to (length, static_cast<Real>(0.));
1252 
1253 #if LIBMESH_DIM == 1
1254  return TypeVector<T>(_coords[0]/length);
1255 #endif
1256 
1257 #if LIBMESH_DIM == 2
1258  return TypeVector<T>(_coords[0]/length,
1259  _coords[1]/length);
1260 #endif
1261 
1262 #if LIBMESH_DIM == 3
1263  return TypeVector<T>(_coords[0]/length,
1264  _coords[1]/length,
1265  _coords[2]/length);
1266 #endif
1267 
1268 }
1269 
1270 template <typename T>
1271 struct CompareTypes<TypeVector<T>, TypeVector<T>>
1272 {
1274 };
1275 
1276 template <typename T, typename T2>
1277 struct CompareTypes<TypeVector<T>, TypeVector<T2>>
1278 {
1280 };
1281 
1282 template <typename T, typename T2, typename std::enable_if<ScalarTraits<T>::value, int>::type = 0>
1284 outer_product(const T & a, const TypeVector<T2> & b)
1285 {
1287  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
1288  ret(i) = a * libmesh_conj(b(i));
1289 
1290  return ret;
1291 }
1292 
1293 template <typename T, typename T2, typename std::enable_if<ScalarTraits<T2>::value, int>::type = 0>
1294 TypeVector<typename CompareTypes<T, T2>::supertype>
1295 outer_product(const TypeVector<T> & a, const T2 & b)
1296 {
1298  const auto conj_b = libmesh_conj(b);
1299  for (unsigned int i = 0; i < LIBMESH_DIM; i++)
1300  ret(i) = a(i) * conj_b;
1301 
1302  return ret;
1303 }
1304 } // namespace libMesh
1305 
1306 namespace std
1307 {
1308 template <typename T>
1309 auto norm(const libMesh::TypeVector<T> & vector) -> decltype(std::norm(T()))
1310 {
1311  return vector.norm();
1312 }
1313 } // namespace std
1314 
1315 #ifdef LIBMESH_HAVE_METAPHYSICL
1316 namespace MetaPhysicL
1317 {
1318 template <typename T>
1319 struct RawType<libMesh::TypeVector<T>>
1320 {
1322 
1324  {
1325  value_type ret;
1326  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1327  ret(i) = raw_value(in(i));
1328 
1329  return ret;
1330  }
1331 };
1332 
1333 template <typename T, typename U>
1334 struct ReplaceAlgebraicType<libMesh::TypeVector<T>, U>
1335 {
1336  typedef U type;
1337 };
1338 }
1339 #endif
1340 
1341 #endif // LIBMESH_TYPE_VECTOR_H
T _coords[LIBMESH_DIM]
TypeVector 的坐标。
Definition: type_vector.h:565
bool operator==(const TypeVector< T > &rhs) const
判断两个向量的每个分量是否相等。
Definition: type_vector.h:1155
const TypeVector< T > & operator*=(const T &)
将该向量乘以标量值。
Definition: type_vector.h:926
const T & slice(const unsigned int i) const
Definition: type_vector.h:262
TypeTensor< typename CompareTypes< T, T2 >::supertype > outer_product(const TypeVector< T > &a, const TypeVector< T2 > &b)
Definition: type_tensor.h:1483
auto norm() const -> decltype(std::norm(T()))
返回向量的模,即元素平方和的平方根。
Definition: type_vector.h:1055
T libmesh_conj(T a)
void subtract_scaled(const TypeVector< T2 > &, const T &)
从该向量中减去另一个向量的缩放值,不创建临时对象。
Definition: type_vector.h:851
CompareTypes< T, T2 >::supertype contract(const TypeVector< T2 > &) const
返回 TypeVector::operator*() 的结果。
Definition: type_vector.h:1022
libMesh::TypeVector< typename RawType< T >::value_type > value_type
Definition: type_vector.h:1321
static constexpr Real TOLERANCE
bool is_zero() const
判断向量的所有值是否为零。
Definition: type_vector.h:1095
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Calls cross_norm_sq() and takes the square root of the result.
Definition: type_vector.h:1238
static value_type value(const libMesh::TypeVector< T > &in)
Definition: type_vector.h:1323
TypeVector< typename CompareTypes< T, T2 >::supertype > supertype
Definition: type_vector.h:1279
unsigned int index_type
辅助 typedef 用于泛型索引编程。
Definition: type_vector.h:211
TypeVector< T > operator-() const
返回该向量的相反向量的副本。
Definition: type_vector.h:861
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的平方根。
Definition: type_vector.h:88
void write_unformatted(std::ostream &out_stream, const bool newline=true) const
无格式输出到流 out。只是将向量的元素用空格分隔打印出来。 默认情况下,还会打印一个换行符,但可以通过 newline 参数来控制这个行为。
Definition: type_vector.C:75
TypeVector< typename CompareTypes< T, T2 >::supertype > operator+(const TypeVector< T2 > &) const
向量相加。
Definition: type_vector.h:716
const TypeVector< T > & operator+=(const TypeVector< T2 > &)
向量相加。
Definition: type_vector.h:741
This class defines a tensor in LIBMESH_DIM dimensional space of type T.
Definition: tensor_tools.h:36
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的绝对值。
Definition: type_vector.h:112
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:74
T cross_norm_sq(const TypeVector< T > &b, const TypeVector< T > &c)
Compute |b x c|^2 without creating the extra temporary produced by calling b.cross(c).norm_sq().
Definition: type_vector.h:1217
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1194
void libmesh_ignore(const Args &...)
auto norm_sq() const -> decltype(std::norm(T()))
返回向量的模的平方,即元素模的平方和。
Definition: type_vector.h:1074
void zero()
将向量的所有分量设置为零。
Definition: type_vector.h:1064
void subtract(const TypeVector< T2 > &)
从另一个向量中减去该向量,不创建临时对象。
Definition: type_vector.h:840
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeNTensor< N, typename CompareTypes< Scalar, T >::supertype > >::type operator*(const Scalar &, const TypeNTensor< N, T > &)
T value_type
辅助 typedef 用于 C++98 泛型编程。
Definition: type_vector.h:206
void print(std::ostream &os=libMesh::out) const
格式化输出,默认输出到 libMesh::out 流。
Definition: type_vector.C:47
bool operator!=(const TypeVector< T > &rhs) const
判断两个向量是否不相等。
Definition: type_vector.h:1177
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
计算该向量与另一个向量的叉积。
Definition: type_vector.h:1032
该类定义了一个在 LIBMESH_DIM 维度空间中类型为 T 的向量。
Definition: tensor_tools.h:34
OStreamProxy out
const TypeVector< T > & operator/=(const T &)
将该向量的每个分量除以标量值。
Definition: type_vector.h:983
TypeVector< T > unit() const
返回指向该向量方向的单位向量。
Definition: type_vector.h:1246
TypeVector()
空构造函数。将向量初始化为 LIBMESH_DIM 维中的零向量。
Definition: type_vector.h:577
const TypeVector< T > & operator-=(const TypeVector< T2 > &)
从该向量中减去另一个向量。
Definition: type_vector.h:828
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
判断两个向量是否在绝对容差 tol 内相等。
Definition: type_vector.h:1105
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator/(const Scalar &) const
将该向量的每个分量除以标量值。
Definition: type_vector.h:954
void assign(const TypeVector< T2 > &)
在不创建临时对象的情况下将向量赋值给该向量。
Definition: type_vector.h:682
const T & operator()(const unsigned int i) const
返回向量中的第 i 个分量的常量引用。
Definition: type_vector.h:692
T & slice(const unsigned int i)
Definition: type_vector.h:271
Eigen::Matrix< MetaPhysicL::DualNumber< T, D, asd >, Eigen::Dynamic, 1 > ADRealEigenVector
使用MetaPhysicL::DualNumber类型的Eigen矩阵来定义自动微分实数向量类型。
Definition: type_vector.h:52
ADRealEigenVector< T, D, asd > norm(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的范数。
Definition: type_vector.h:64
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector & >::type operator=(const Scalar &libmesh_dbg_var(p))
赋值操作符,用于将标量值赋给该向量以清零。
Definition: type_vector.h:252
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
判断两个向量是否在相对容差 tol 内相等。
Definition: type_vector.h:1130
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeVector< typename CompareTypes< T, Scalar >::supertype > >::type operator*(const Scalar &) const
将该向量乘以标量值。
Definition: type_vector.h:889
~TypeVector()=default
析构函数。