type_tensor.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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_TENSOR_H
21 #define LIBMESH_TYPE_TENSOR_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/type_vector.h"
26 
27 // C++ includes
28 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
29 #include <cmath>
30 #include <tuple>
31 
32 namespace libMesh
33 {
34 
35 // Forward declarations
36 template <typename T> class TypeTensorColumn;
37 template <typename T> class ConstTypeTensorColumn;
38 template <unsigned int N, typename T> class TypeNTensor;
39 
47 template <typename T>
48 class TypeTensor
49 {
50  template <typename T2>
51  friend class TypeTensor;
52 
53 public:
58  TypeTensor ();
59 
60 protected:
67  explicit TypeTensor (const T xx,
68  const T xy=0,
69  const T xz=0,
70  const T yx=0,
71  const T yy=0,
72  const T yz=0,
73  const T zx=0,
74  const T zy=0,
75  const T zz=0);
76 
77 
81  template <typename Scalar>
82  explicit TypeTensor (const Scalar xx,
83  const Scalar xy=0,
84  const Scalar xz=0,
85  const Scalar yx=0,
86  const Scalar yy=0,
87  const Scalar yz=0,
88  const Scalar zx=0,
89  const Scalar zy=0,
90  typename
92  const Scalar>::type zz=0);
93 
99  template<typename T2>
100  TypeTensor(const TypeVector<T2> & vx);
101 
102  template<typename T2>
103  TypeTensor(const TypeVector<T2> & vx,
104  const TypeVector<T2> & vy);
105 
106  template<typename T2>
107  TypeTensor(const TypeVector<T2> & vx,
108  const TypeVector<T2> & vy,
109  const TypeVector<T2> & vz);
110 
111 public:
112 
116  typedef T value_type;
117 
121  typedef std::tuple<unsigned int, unsigned int> index_type;
122 
126  template<typename T2>
127  TypeTensor(const TypeTensor<T2> & p);
128 
132  ~TypeTensor();
133 
137  template<typename T2>
138  void assign (const TypeTensor<T2> &);
139 
145  template <typename Scalar>
146  typename boostcopy::enable_if_c<
148  TypeTensor &>::type
149  operator = (const Scalar & libmesh_dbg_var(p))
150  { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; }
151 
155  const T & operator () (const unsigned int i, const unsigned int j) const;
156 
160  T & operator () (const unsigned int i, const unsigned int j);
161 
165  ConstTypeTensorColumn<T> slice (const unsigned int i) const;
166 
170  TypeTensorColumn<T> slice (const unsigned int i);
171 
175  TypeVector<T> row(const unsigned int r) const;
176 
182  template<typename T2>
184  operator + (const TypeTensor<T2> &) const;
185 
191  template<typename T2>
192  const TypeTensor<T> & operator += (const TypeTensor<T2> &);
193 
197  template<typename T2>
198  void add (const TypeTensor<T2> &);
199 
203  template <typename T2>
204  void add_scaled (const TypeTensor<T2> &, const T);
205 
211  template<typename T2>
213  operator - (const TypeTensor<T2> &) const;
214 
220  template<typename T2>
221  const TypeTensor<T> & operator -= (const TypeTensor<T2> &);
222 
226  template<typename T2>
227  void subtract (const TypeTensor<T2> &);
228 
233  template <typename T2>
234  void subtract_scaled (const TypeTensor<T2> &, const T);
235 
239  TypeTensor<T> operator - () const;
240 
246  template <typename Scalar>
247  auto
248  operator * (const Scalar scalar) const -> typename boostcopy::enable_if_c<
251 
257  template <typename Scalar, typename boostcopy::enable_if_c<
258  ScalarTraits<Scalar>::value, int>::type = 0>
259  const TypeTensor<T> & operator *= (const Scalar factor)
260  {
261  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
262  _coords[i] *= factor;
263 
264  return *this;
265  }
266 
272  template <typename Scalar>
273  typename boostcopy::enable_if_c<
276  operator / (const Scalar) const;
277 
283  const TypeTensor<T> & operator /= (const T);
284 
291  template <typename T2>
293 
299  template <typename T2>
300  const TypeTensor<T> & operator *= (const TypeTensor<T2> &);
301 
311  template <typename T2>
313  contract (const TypeTensor<T2> &) const;
314 
321  template <typename T2>
323  operator * (const TypeVector<T2> &) const;
324 
328  TypeTensor<T> transpose() const;
329 
333  TypeTensor<T> inverse() const;
334 
341  void solve(const TypeVector<T> & b, TypeVector<T> & x) const;
342 
349 #ifdef LIBMESH_ENABLE_DEPRECATED
350  Real size() const;
351 #endif
352 
357  Real norm() const;
358 
365 #ifdef LIBMESH_ENABLE_DEPRECATED
366  Real size_sq() const;
367 #endif
368 
373  Real norm_sq() const;
374 
381  T det() const;
382 
386  T tr() const;
387 
391  void zero();
392 
396  bool operator == (const TypeTensor<T> & rhs) const;
397 
402  bool operator < (const TypeTensor<T> & rhs) const;
403 
407  bool operator > (const TypeTensor<T> & rhs) const;
408 
412  void print(std::ostream & os = libMesh::out) const;
413 
421  friend std::ostream & operator << (std::ostream & os, const TypeTensor<T> & t)
422  {
423  t.print(os);
424  return os;
425  }
426 
431  void write_unformatted (std::ostream & out, const bool newline = true) const;
432 
433 protected:
434 
438  T _coords[LIBMESH_DIM*LIBMESH_DIM];
439 };
440 
441 
442 template <typename T>
443 class TypeTensorColumn
444 {
445 public:
446 
448  unsigned int j) :
449  _tensor(&tensor), _j(j) {}
450 
455  T & operator () (const unsigned int i)
456  { return (*_tensor)(i,_j); }
457 
458  T & slice (const unsigned int i)
459  { return (*_tensor)(i,_j); }
460 
465  {
466  for (unsigned int i=0; i != LIBMESH_DIM; ++i)
467  (*this)(i) = rhs(i);
468  return *this;
469  }
470 
471 private:
473  const unsigned int _j;
474 };
475 
476 
477 template <typename T>
479 {
480 public:
481 
483  unsigned int j) :
484  _tensor(&tensor), _j(j) {}
485 
489  const T & operator () (const unsigned int i) const
490  { return (*_tensor)(i,_j); }
491 
492  const T & slice (const unsigned int i) const
493  { return (*_tensor)(i,_j); }
494 
495 private:
497  const unsigned int _j;
498 };
499 
500 //------------------------------------------------------
501 // Inline functions
502 template <typename T>
503 inline
505 {
506  _coords[0] = 0;
507 
508 #if LIBMESH_DIM > 1
509  _coords[1] = 0;
510  _coords[2] = 0;
511  _coords[3] = 0;
512 #endif
513 
514 #if LIBMESH_DIM > 2
515  _coords[4] = 0;
516  _coords[5] = 0;
517  _coords[6] = 0;
518  _coords[7] = 0;
519  _coords[8] = 0;
520 #endif
521 }
522 
523 
524 
525 template <typename T>
526 inline
528  const T xy,
529  const T xz,
530  const T yx,
531  const T yy,
532  const T yz,
533  const T zx,
534  const T zy,
535  const T zz)
536 {
537  _coords[0] = xx;
538 
539 #if LIBMESH_DIM == 2
540  _coords[1] = xy;
541  _coords[2] = yx;
542  _coords[3] = yy;
543 #elif LIBMESH_DIM == 1
544  libmesh_assert_equal_to (xy, 0);
545  libmesh_assert_equal_to (yx, 0);
546  libmesh_assert_equal_to (yy, 0);
547 #endif
548 
549 #if LIBMESH_DIM == 3
550  _coords[1] = xy;
551  _coords[2] = xz;
552  _coords[3] = yx;
553  _coords[4] = yy;
554  _coords[5] = yz;
555  _coords[6] = zx;
556  _coords[7] = zy;
557  _coords[8] = zz;
558 #else
559  libmesh_assert_equal_to (xz, 0);
560  libmesh_assert_equal_to (yz, 0);
561  libmesh_assert_equal_to (zx, 0);
562  libmesh_assert_equal_to (zy, 0);
563  libmesh_assert_equal_to (zz, 0);
564 #endif
565 }
566 
567 
568 template <typename T>
569 template <typename Scalar>
570 inline
571 TypeTensor<T>::TypeTensor (const Scalar xx,
572  const Scalar xy,
573  const Scalar xz,
574  const Scalar yx,
575  const Scalar yy,
576  const Scalar yz,
577  const Scalar zx,
578  const Scalar zy,
579  typename
581  const Scalar>::type zz)
582 {
583  _coords[0] = xx;
584 
585 #if LIBMESH_DIM == 2
586  _coords[1] = xy;
587  _coords[2] = yx;
588  _coords[3] = yy;
589 #endif
590 
591 #if LIBMESH_DIM == 3
592  _coords[1] = xy;
593  _coords[2] = xz;
594  _coords[3] = yx;
595  _coords[4] = yy;
596  _coords[5] = yz;
597  _coords[6] = zx;
598  _coords[7] = zy;
599  _coords[8] = zz;
600 #endif
601 }
602 
603 
604 
605 template <typename T>
606 template<typename T2>
607 inline
609 {
610  // copy the nodes from vector p to me
611  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
612  _coords[i] = p._coords[i];
613 }
614 
615 
616 template <typename T>
617 template <typename T2>
619 {
620  libmesh_assert_equal_to (LIBMESH_DIM, 1);
621  _coords[0] = vx(0);
622 }
623 
624 template <typename T>
625 template <typename T2>
627  const TypeVector<T2> & vy)
628 {
629  libmesh_assert_equal_to (LIBMESH_DIM, 2);
630  _coords[0] = vx(0);
631  _coords[1] = vx(1);
632  _coords[2] = vy(0);
633  _coords[3] = vy(1);
634 }
635 
636 template <typename T>
637 template <typename T2>
639  const TypeVector<T2> & vy,
640  const TypeVector<T2> & vz)
641 {
642  libmesh_assert_equal_to (LIBMESH_DIM, 3);
643  _coords[0] = vx(0);
644  _coords[1] = vx(1);
645  _coords[2] = vx(2);
646  _coords[3] = vy(0);
647  _coords[4] = vy(1);
648  _coords[5] = vy(2);
649  _coords[6] = vz(0);
650  _coords[7] = vz(1);
651  _coords[8] = vz(2);
652 }
653 
654 
655 
656 
657 template <typename T>
658 inline
660 {
661 }
662 
663 
664 
665 template <typename T>
666 template<typename T2>
667 inline
669 {
670  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
671  _coords[i] = p._coords[i];
672 }
673 
674 
675 
676 template <typename T>
677 inline
678 const T & TypeTensor<T>::operator () (const unsigned int i,
679  const unsigned int j) const
680 {
681  libmesh_assert_less (i, 3);
682  libmesh_assert_less (j, 3);
683 
684 #if LIBMESH_DIM < 3
685  const static T my_zero = 0;
686  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
687  return my_zero;
688 #endif
689 
690  return _coords[i*LIBMESH_DIM+j];
691 }
692 
693 
694 
695 template <typename T>
696 inline
697 T & TypeTensor<T>::operator () (const unsigned int i,
698  const unsigned int j)
699 {
700 #if LIBMESH_DIM < 3
701 
702  if (i >= LIBMESH_DIM || j >= LIBMESH_DIM)
703  libmesh_error_msg("ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!");
704 
705 #endif
706 
707  libmesh_assert_less (i, LIBMESH_DIM);
708  libmesh_assert_less (j, LIBMESH_DIM);
709 
710  return _coords[i*LIBMESH_DIM+j];
711 }
712 
713 
714 template <typename T>
715 inline
717 TypeTensor<T>::slice (const unsigned int i) const
718 {
719  libmesh_assert_less (i, LIBMESH_DIM);
720  return ConstTypeTensorColumn<T>(*this, i);
721 }
722 
723 
724 template <typename T>
725 inline
727 TypeTensor<T>::slice (const unsigned int i)
728 {
729  libmesh_assert_less (i, LIBMESH_DIM);
730  return TypeTensorColumn<T>(*this, i);
731 }
732 
733 
734 template <typename T>
735 inline
737 TypeTensor<T>::row(const unsigned int r) const
738 {
739  TypeVector<T> return_vector;
740 
741  for (unsigned int j=0; j<LIBMESH_DIM; j++)
742  return_vector._coords[j] = _coords[r*LIBMESH_DIM + j];
743 
744  return return_vector;
745 }
746 
747 template <typename T>
748 template<typename T2>
749 inline
752 {
753 
754 #if LIBMESH_DIM == 1
755  return TypeTensor(_coords[0] + p._coords[0]);
756 #endif
757 
758 #if LIBMESH_DIM == 2
759  return TypeTensor(_coords[0] + p._coords[0],
760  _coords[1] + p._coords[1],
761  0.,
762  _coords[2] + p._coords[2],
763  _coords[3] + p._coords[3]);
764 #endif
765 
766 #if LIBMESH_DIM == 3
767  return TypeTensor(_coords[0] + p._coords[0],
768  _coords[1] + p._coords[1],
769  _coords[2] + p._coords[2],
770  _coords[3] + p._coords[3],
771  _coords[4] + p._coords[4],
772  _coords[5] + p._coords[5],
773  _coords[6] + p._coords[6],
774  _coords[7] + p._coords[7],
775  _coords[8] + p._coords[8]);
776 #endif
777 
778 }
779 
780 
781 
782 template <typename T>
783 template<typename T2>
784 inline
786 {
787  this->add (p);
788 
789  return *this;
790 }
791 
792 
793 
794 template <typename T>
795 template<typename T2>
796 inline
798 {
799  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
800  _coords[i] += p._coords[i];
801 }
802 
803 
804 
805 template <typename T>
806 template <typename T2>
807 inline
808 void TypeTensor<T>::add_scaled (const TypeTensor<T2> & p, const T factor)
809 {
810  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
811  _coords[i] += factor*p._coords[i];
812 
813 }
814 
815 
816 
817 template <typename T>
818 template<typename T2>
819 inline
822 {
823 
824 #if LIBMESH_DIM == 1
825  return TypeTensor(_coords[0] - p._coords[0]);
826 #endif
827 
828 #if LIBMESH_DIM == 2
829  return TypeTensor(_coords[0] - p._coords[0],
830  _coords[1] - p._coords[1],
831  0.,
832  _coords[2] - p._coords[2],
833  _coords[3] - p._coords[3]);
834 #endif
835 
836 #if LIBMESH_DIM == 3
837  return TypeTensor(_coords[0] - p._coords[0],
838  _coords[1] - p._coords[1],
839  _coords[2] - p._coords[2],
840  _coords[3] - p._coords[3],
841  _coords[4] - p._coords[4],
842  _coords[5] - p._coords[5],
843  _coords[6] - p._coords[6],
844  _coords[7] - p._coords[7],
845  _coords[8] - p._coords[8]);
846 #endif
847 
848 }
849 
850 
851 
852 template <typename T>
853 template <typename T2>
854 inline
856 {
857  this->subtract (p);
858 
859  return *this;
860 }
861 
862 
863 
864 template <typename T>
865 template <typename T2>
866 inline
868 {
869  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
870  _coords[i] -= p._coords[i];
871 }
872 
873 
874 
875 template <typename T>
876 template <typename T2>
877 inline
878 void TypeTensor<T>::subtract_scaled (const TypeTensor<T2> & p, const T factor)
879 {
880  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
881  _coords[i] -= factor*p._coords[i];
882 }
883 
884 
885 
886 template <typename T>
887 inline
889 {
890 
891 #if LIBMESH_DIM == 1
892  return TypeTensor(-_coords[0]);
893 #endif
894 
895 #if LIBMESH_DIM == 2
896  return TypeTensor(-_coords[0],
897  -_coords[1],
898  -_coords[2],
899  -_coords[3]);
900 #endif
901 
902 #if LIBMESH_DIM == 3
903  return TypeTensor(-_coords[0],
904  -_coords[1],
905  -_coords[2],
906  -_coords[3],
907  -_coords[4],
908  -_coords[5],
909  -_coords[6],
910  -_coords[7],
911  -_coords[8]);
912 #endif
913 
914 }
915 
916 
917 
918 template <typename T>
919 template <typename Scalar>
920 inline
921 auto
922 TypeTensor<T>::operator * (const Scalar factor) const -> typename boostcopy::enable_if_c<
925 {
926  typedef decltype((*this)(0, 0) * factor) TS;
927 
928 
929 #if LIBMESH_DIM == 1
930  return TypeTensor<TS>(_coords[0]*factor);
931 #endif
932 
933 #if LIBMESH_DIM == 2
934  return TypeTensor<TS>(_coords[0]*factor,
935  _coords[1]*factor,
936  _coords[2]*factor,
937  _coords[3]*factor);
938 #endif
939 
940 #if LIBMESH_DIM == 3
941  return TypeTensor<TS>(_coords[0]*factor,
942  _coords[1]*factor,
943  _coords[2]*factor,
944  _coords[3]*factor,
945  _coords[4]*factor,
946  _coords[5]*factor,
947  _coords[6]*factor,
948  _coords[7]*factor,
949  _coords[8]*factor);
950 #endif
951 }
952 
953 
954 
955 template <typename T, typename Scalar>
956 inline
957 typename boostcopy::enable_if_c<
960 operator * (const Scalar factor,
961  const TypeTensor<T> & t)
962 {
963  return t * factor;
964 }
965 
966 template <typename T>
967 template <typename Scalar>
968 inline
969 typename boostcopy::enable_if_c<
971  TypeTensor<typename CompareTypes<T, Scalar>::supertype>>::type
972 TypeTensor<T>::operator / (const Scalar factor) const
973 {
974  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
975 
976  typedef typename CompareTypes<T, Scalar>::supertype TS;
977 
978 #if LIBMESH_DIM == 1
979  return TypeTensor<TS>(_coords[0]/factor);
980 #endif
981 
982 #if LIBMESH_DIM == 2
983  return TypeTensor<TS>(_coords[0]/factor,
984  _coords[1]/factor,
985  _coords[2]/factor,
986  _coords[3]/factor);
987 #endif
988 
989 #if LIBMESH_DIM == 3
990  return TypeTensor<TS>(_coords[0]/factor,
991  _coords[1]/factor,
992  _coords[2]/factor,
993  _coords[3]/factor,
994  _coords[4]/factor,
995  _coords[5]/factor,
996  _coords[6]/factor,
997  _coords[7]/factor,
998  _coords[8]/factor);
999 #endif
1000 
1001 }
1002 
1003 
1004 
1005 template <typename T>
1006 inline
1008 {
1009 #if LIBMESH_DIM == 1
1010  return TypeTensor(_coords[0]);
1011 #endif
1012 
1013 #if LIBMESH_DIM == 2
1014  return TypeTensor(_coords[0],
1015  _coords[2],
1016  _coords[1],
1017  _coords[3]);
1018 #endif
1019 
1020 #if LIBMESH_DIM == 3
1021  return TypeTensor(_coords[0],
1022  _coords[3],
1023  _coords[6],
1024  _coords[1],
1025  _coords[4],
1026  _coords[7],
1027  _coords[2],
1028  _coords[5],
1029  _coords[8]);
1030 #endif
1031 }
1032 
1033 
1034 
1035 template <typename T>
1036 inline
1038 {
1039 #if LIBMESH_DIM == 1
1040  if (_coords[0] == static_cast<T>(0.))
1041  libmesh_convergence_failure();
1042  return TypeTensor(1. / _coords[0]);
1043 #endif
1044 
1045 #if LIBMESH_DIM == 2
1046  // Get convenient reference to this.
1047  const TypeTensor<T> & A = *this;
1048 
1049  // Use temporary variables, avoid multiple accesses
1050  T a = A(0,0), b = A(0,1),
1051  c = A(1,0), d = A(1,1);
1052 
1053  // Make sure det = ad - bc is not zero
1054  T my_det = a*d - b*c;
1055 
1056  if (my_det == static_cast<T>(0.))
1057  libmesh_convergence_failure();
1058 
1059  return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det);
1060 #endif
1061 
1062 #if LIBMESH_DIM == 3
1063  // Get convenient reference to this.
1064  const TypeTensor<T> & A = *this;
1065 
1066  T a11 = A(0,0), a12 = A(0,1), a13 = A(0,2),
1067  a21 = A(1,0), a22 = A(1,1), a23 = A(1,2),
1068  a31 = A(2,0), a32 = A(2,1), a33 = A(2,2);
1069 
1070  T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13);
1071 
1072  if (my_det == static_cast<T>(0.))
1073  libmesh_convergence_failure();
1074 
1075  // Inline comment characters are for lining up columns.
1076  return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det,
1077  -(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det,
1078  (a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det);
1079 #endif
1080 }
1081 
1082 
1083 
1084 template <typename T>
1085 inline
1087 {
1088 #if LIBMESH_DIM == 1
1089  if (_coords[0] == static_cast<T>(0.))
1090  libmesh_convergence_failure();
1091  x(0) = b(0) / _coords[0];
1092 #endif
1093 
1094 #if LIBMESH_DIM == 2
1095  T my_det = _coords[0]*_coords[3] - _coords[1]*_coords[2];
1096 
1097  if (my_det == static_cast<T>(0.))
1098  libmesh_convergence_failure();
1099 
1100  T my_det_inv = 1./my_det;
1101 
1102  x(0) = my_det_inv*( _coords[3]*b(0) - _coords[1]*b(1));
1103  x(1) = my_det_inv*(-_coords[2]*b(0) + _coords[0]*b(1));
1104 #endif
1105 
1106 #if LIBMESH_DIM == 3
1107  T my_det =
1108  // a11*(a33 *a22 - a32 *a23)
1109  _coords[0]*(_coords[8]*_coords[4] - _coords[7]*_coords[5])
1110  // -a21*(a33 *a12 - a32 *a13)
1111  -_coords[3]*(_coords[8]*_coords[1] - _coords[7]*_coords[2]) +
1112  // +a31*(a23 *a12 - a22 *a13)
1113  +_coords[6]*(_coords[5]*_coords[1] - _coords[4]*_coords[2]);
1114 
1115  if (my_det == static_cast<T>(0.))
1116  libmesh_convergence_failure();
1117 
1118  T my_det_inv = 1./my_det;
1119  x(0) = my_det_inv*((_coords[8]*_coords[4] - _coords[7]*_coords[5])*b(0) -
1120  (_coords[8]*_coords[1] - _coords[7]*_coords[2])*b(1) +
1121  (_coords[5]*_coords[1] - _coords[4]*_coords[2])*b(2));
1122 
1123  x(1) = my_det_inv*((_coords[6]*_coords[5] - _coords[8]*_coords[3])*b(0) +
1124  (_coords[8]*_coords[0] - _coords[6]*_coords[2])*b(1) -
1125  (_coords[5]*_coords[0] - _coords[3]*_coords[2])*b(2));
1126 
1127  x(2) = my_det_inv*((_coords[7]*_coords[3] - _coords[6]*_coords[4])*b(0) -
1128  (_coords[7]*_coords[0] - _coords[6]*_coords[1])*b(1) +
1129  (_coords[4]*_coords[0] - _coords[3]*_coords[1])*b(2));
1130 #endif
1131 }
1132 
1133 
1134 
1135 template <typename T>
1136 inline
1138 {
1139  libmesh_assert_not_equal_to (factor, static_cast<T>(0.));
1140 
1141  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1142  _coords[i] /= factor;
1143 
1144  return *this;
1145 }
1146 
1147 
1148 
1149 
1150 template <typename T>
1151 template <typename T2>
1152 inline
1155 {
1157  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1158  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1159  returnval(i) += (*this)(i,j)*p(j);
1160 
1161  return returnval;
1162 }
1163 
1164 
1165 
1166 template <typename T>
1167 template <typename T2>
1168 inline
1170 {
1171  TypeTensor<T> returnval;
1172  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1173  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1174  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1175  returnval(i,j) += (*this)(i,k)*p(k,j);
1176 
1177  return returnval;
1178 }
1179 
1180 template <typename T>
1181 template <typename T2>
1182 inline
1184 {
1185  TypeTensor<T> temp;
1186  for (unsigned int i=0; i<LIBMESH_DIM; i++)
1187  for (unsigned int j=0; j<LIBMESH_DIM; j++)
1188  for (unsigned int k=0; k<LIBMESH_DIM; k++)
1189  temp(i,j) += (*this)(i,k)*p(k,j);
1190 
1191  this->assign(temp);
1192  return *this;
1193 }
1194 
1195 
1200 template <typename T>
1201 template <typename T2>
1202 inline
1205 {
1206  typename CompareTypes<T,T2>::supertype sum = 0.;
1207  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1208  sum += _coords[i]*t._coords[i];
1209  return sum;
1210 }
1211 
1212 
1213 
1214 #ifdef LIBMESH_ENABLE_DEPRECATED
1215 template <typename T>
1216 inline
1218 {
1219  libmesh_deprecated();
1220  return this->norm();
1221 }
1222 #endif
1223 
1224 
1225 
1226 template <typename T>
1227 inline
1229 {
1230  return std::sqrt(this->norm_sq());
1231 }
1232 
1233 
1234 
1235 template <typename T>
1236 inline
1238 {
1239 #if LIBMESH_DIM == 1
1240  return _coords[0];
1241 #endif
1242 
1243 #if LIBMESH_DIM == 2
1244  return (_coords[0] * _coords[3]
1245  - _coords[1] * _coords[2]);
1246 #endif
1247 
1248 #if LIBMESH_DIM == 3
1249  return (_coords[0] * _coords[4] * _coords[8]
1250  + _coords[1] * _coords[5] * _coords[6]
1251  + _coords[2] * _coords[3] * _coords[7]
1252  - _coords[0] * _coords[5] * _coords[7]
1253  - _coords[1] * _coords[3] * _coords[8]
1254  - _coords[2] * _coords[4] * _coords[6]);
1255 #endif
1256 }
1257 
1258 template <typename T>
1259 inline
1261 {
1262 #if LIBMESH_DIM == 1
1263  return _coords[0];
1264 #endif
1265 
1266 #if LIBMESH_DIM == 2
1267  return _coords[0] + _coords[3];
1268 #endif
1269 
1270 #if LIBMESH_DIM == 3
1271  return _coords[0] + _coords[4] + _coords[8];
1272 #endif
1273 }
1274 
1275 template <typename T>
1276 inline
1278 {
1279  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1280  _coords[i] = 0.;
1281 }
1282 
1283 
1284 
1285 #ifdef LIBMESH_ENABLE_DEPRECATED
1286 template <typename T>
1287 inline
1289 {
1290  libmesh_deprecated();
1291  return this->norm_sq();
1292 }
1293 #endif
1294 
1295 
1296 
1297 template <typename T>
1298 inline
1300 {
1301  Real sum = 0.;
1302  for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++)
1303  sum += TensorTools::norm_sq(_coords[i]);
1304  return sum;
1305 }
1306 
1307 
1308 
1309 template <typename T>
1310 inline
1312 {
1313 #if LIBMESH_DIM == 1
1314  return (std::abs(_coords[0] - rhs._coords[0])
1315  < TOLERANCE);
1316 #endif
1317 
1318 #if LIBMESH_DIM == 2
1319  return ((std::abs(_coords[0] - rhs._coords[0]) +
1320  std::abs(_coords[1] - rhs._coords[1]) +
1321  std::abs(_coords[2] - rhs._coords[2]) +
1322  std::abs(_coords[3] - rhs._coords[3]))
1323  < 4.*TOLERANCE);
1324 #endif
1325 
1326 #if LIBMESH_DIM == 3
1327  return ((std::abs(_coords[0] - rhs._coords[0]) +
1328  std::abs(_coords[1] - rhs._coords[1]) +
1329  std::abs(_coords[2] - rhs._coords[2]) +
1330  std::abs(_coords[3] - rhs._coords[3]) +
1331  std::abs(_coords[4] - rhs._coords[4]) +
1332  std::abs(_coords[5] - rhs._coords[5]) +
1333  std::abs(_coords[6] - rhs._coords[6]) +
1334  std::abs(_coords[7] - rhs._coords[7]) +
1335  std::abs(_coords[8] - rhs._coords[8]))
1336  < 9.*TOLERANCE);
1337 #endif
1338 
1339 }
1340 
1341 } // namespace libMesh
1342 
1343 #endif // LIBMESH_TYPE_TENSOR_H
void write_unformatted(std::ostream &out, const bool newline=true) const
Definition: type_tensor.C:83
const T & operator()(const unsigned int i) const
Definition: type_tensor.h:489
TypeTensorColumn< T > & operator=(const TypeVector< T > &rhs)
Definition: type_tensor.h:464
T _coords[LIBMESH_DIM]
Definition: type_vector.h:409
double abs(double a)
const unsigned int _j
Definition: type_tensor.h:473
const T & slice(const unsigned int i) const
Definition: type_tensor.h:492
TypeTensor< T > * _tensor
Definition: type_tensor.h:472
T & operator()(const unsigned int i)
Definition: type_tensor.h:455
void add_scaled(const TypeTensor< T2 > &, const T)
Definition: type_tensor.h:808
std::tuple< unsigned int, unsigned int > index_type
Definition: type_tensor.h:121
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< typename CompareTypes< T, Scalar >::supertype > >::type operator/(const Scalar) const
Definition: type_tensor.h:972
TypeTensor< T > transpose() const
Definition: type_tensor.h:1007
ConstTypeTensorColumn(const TypeTensor< T > &tensor, unsigned int j)
Definition: type_tensor.h:482
void solve(const TypeVector< T > &b, TypeVector< T > &x) const
Definition: type_tensor.h:1086
TypeVector< T > row(const unsigned int r) const
Definition: type_tensor.h:737
static const Real TOLERANCE
T _coords[LIBMESH_DIM *LIBMESH_DIM]
Definition: type_tensor.h:438
const TypeTensor< T > * _tensor
Definition: type_tensor.h:496
TypeTensor< typename CompareTypes< T, T2 >::supertype > operator+(const TypeTensor< T2 > &) const
Definition: type_tensor.h:751
TypeTensorColumn(TypeTensor< T > &tensor, unsigned int j)
Definition: type_tensor.h:447
auto operator*(const Scalar scalar) const -> typename boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor< decltype(T() *scalar)>>::type
Definition: type_tensor.h:922
TypeTensor< T > operator-() const
Definition: type_tensor.h:888
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeNTensor< N, typename CompareTypes< Scalar, T >::supertype > >::type operator*(const Scalar &, const TypeNTensor< N, T > &)
const TypeTensor< T > & operator/=(const T)
Definition: type_tensor.h:1137
void subtract(const TypeTensor< T2 > &)
Definition: type_tensor.h:867
void assign(const TypeTensor< T2 > &)
Definition: type_tensor.h:668
Real norm() const
Definition: type_tensor.h:1228
void add(const TypeTensor< T2 > &)
Definition: type_tensor.h:797
CompareTypes< T, T2 >::supertype contract(const TypeTensor< T2 > &) const
Definition: type_tensor.h:1204
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const TypeTensor< T > & operator*=(const Scalar factor)
Definition: type_tensor.h:259
TypeTensor< T > inverse() const
Definition: type_tensor.h:1037
static PetscErrorCode Mat * A
boostcopy::enable_if_c< ScalarTraits< Scalar >::value, TypeTensor & >::type operator=(const Scalar &libmesh_dbg_var(p))
Definition: type_tensor.h:149
const TypeTensor< T > & operator-=(const TypeTensor< T2 > &)
Definition: type_tensor.h:855
ConstTypeTensorColumn< T > slice(const unsigned int i) const
Definition: type_tensor.h:717
const TypeTensor< T > & operator+=(const TypeTensor< T2 > &)
Definition: type_tensor.h:785
T & slice(const unsigned int i)
Definition: type_tensor.h:458
bool operator==(const TypeTensor< T > &rhs) const
Definition: type_tensor.h:1311
bool operator>(const TypeTensor< T > &rhs) const
static const bool value
Definition: compare_types.h:58
OStreamProxy out(std::cout)
void print(std::ostream &os=libMesh::out) const
Definition: type_tensor.C:39
const T & operator()(const unsigned int i, const unsigned int j) const
Definition: type_tensor.h:678
Real norm_sq() const
Definition: type_tensor.h:1299
Real size() const
Definition: type_tensor.h:1217
void subtract_scaled(const TypeTensor< T2 > &, const T)
Definition: type_tensor.h:878
Real size_sq() const
Definition: type_tensor.h:1288