$extrastylesheet
00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 #ifndef LIBMESH_TYPE_TENSOR_H 00021 #define LIBMESH_TYPE_TENSOR_H 00022 00023 // Local includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/type_vector.h" 00026 00027 // C++ includes 00028 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC 00029 #include <cmath> 00030 00031 namespace libMesh 00032 { 00033 00034 // Forward declarations 00035 template <typename T> class TypeTensorColumn; 00036 template <typename T> class ConstTypeTensorColumn; 00037 template <unsigned int N, typename T> class TypeNTensor; 00038 00039 00040 00050 template <typename T> 00051 class TypeTensor 00052 { 00053 template <typename T2> 00054 friend class TypeTensor; 00055 00056 protected: 00057 00062 TypeTensor (); 00063 00070 explicit TypeTensor (const T xx, 00071 const T xy=0, 00072 const T xz=0, 00073 const T yx=0, 00074 const T yy=0, 00075 const T yz=0, 00076 const T zx=0, 00077 const T zy=0, 00078 const T zz=0); 00079 00080 00084 template <typename Scalar> 00085 explicit TypeTensor (const Scalar xx, 00086 const Scalar xy=0, 00087 const Scalar xz=0, 00088 const Scalar yx=0, 00089 const Scalar yy=0, 00090 const Scalar yz=0, 00091 const Scalar zx=0, 00092 const Scalar zy=0, 00093 typename 00094 boostcopy::enable_if_c<ScalarTraits<Scalar>::value, 00095 const Scalar>::type zz=0); 00096 00102 template<typename T2> 00103 TypeTensor(const TypeVector<T2>& vx); 00104 00105 template<typename T2> 00106 TypeTensor(const TypeVector<T2>& vx, const TypeVector<T2> &vy); 00107 00108 template<typename T2> 00109 TypeTensor(const TypeVector<T2>& vx, const TypeVector<T2> &vy, const TypeVector<T2> &vz); 00110 00111 public: 00112 00116 template<typename T2> 00117 TypeTensor(const TypeTensor<T2>& p); 00118 00122 ~TypeTensor(); 00123 00127 template<typename T2> 00128 void assign (const TypeTensor<T2> &); 00129 00133 template <typename Scalar> 00134 typename boostcopy::enable_if_c< 00135 ScalarTraits<Scalar>::value, 00136 TypeTensor&>::type 00137 operator = (const Scalar& p) 00138 { libmesh_assert_equal_to (p, Scalar(0)); this->zero(); return *this; } 00139 00143 const T & operator () (const unsigned int i, const unsigned int j) const; 00144 00149 T & operator () (const unsigned int i, const unsigned int j); 00150 00154 ConstTypeTensorColumn<T> slice (const unsigned int i) const; 00155 00159 TypeTensorColumn<T> slice (const unsigned int i); 00160 00164 TypeVector<T> row(const unsigned int r); 00165 00169 template<typename T2> 00170 TypeTensor<typename CompareTypes<T, T2>::supertype> 00171 operator + (const TypeTensor<T2> &) const; 00172 00176 template<typename T2> 00177 const TypeTensor<T> & operator += (const TypeTensor<T2> &); 00178 00182 template<typename T2> 00183 void add (const TypeTensor<T2> &); 00184 00189 template <typename T2> 00190 void add_scaled (const TypeTensor<T2> &, const T); 00191 00195 template<typename T2> 00196 TypeTensor<typename CompareTypes<T, T2>::supertype> 00197 operator - (const TypeTensor<T2> &) const; 00198 00202 template<typename T2> 00203 const TypeTensor<T> & operator -= (const TypeTensor<T2> &); 00204 00208 template<typename T2> 00209 void subtract (const TypeTensor<T2> &); 00210 00215 template <typename T2> 00216 void subtract_scaled (const TypeTensor<T2> &, const T); 00217 00221 TypeTensor<T> operator - () const; 00222 00226 template <typename Scalar> 00227 typename boostcopy::enable_if_c< 00228 ScalarTraits<Scalar>::value, 00229 TypeTensor<typename CompareTypes<T, Scalar>::supertype> >::type 00230 operator * (const Scalar) const; 00231 00235 template <typename Scalar> 00236 const TypeTensor<T> & operator *= (const Scalar); 00237 00241 template <typename Scalar> 00242 typename boostcopy::enable_if_c< 00243 ScalarTraits<Scalar>::value, 00244 TypeTensor<typename CompareTypes<T, Scalar>::supertype> >::type 00245 operator / (const Scalar) const; 00246 00250 const TypeTensor<T> & operator /= (const T); 00251 00256 template <typename T2> 00257 TypeTensor<T> operator * (const TypeTensor<T2> &) const; 00258 00264 template <typename T2> 00265 typename CompareTypes<T,T2>::supertype 00266 contract (const TypeTensor<T2> &) const; 00267 00272 template <typename T2> 00273 TypeVector<typename CompareTypes<T,T2>::supertype> 00274 operator * (const TypeVector<T2> &) const; 00275 00279 TypeTensor<T> transpose() const; 00280 00284 TypeTensor<T> inverse() const; 00285 00290 Real size() const; 00291 00296 Real size_sq() const; 00297 00303 T det() const; 00304 00308 T tr() const; 00309 00313 void zero(); 00314 00318 bool operator == (const TypeTensor<T>& rhs) const; 00319 00324 bool operator < (const TypeTensor<T>& rhs) const; 00325 00330 bool operator > (const TypeTensor<T>& rhs) const; 00331 00335 void print(std::ostream& os = libMesh::out) const; 00336 00341 friend std::ostream& operator << (std::ostream& os, const TypeTensor<T>& t) 00342 { 00343 t.print(os); 00344 return os; 00345 } 00346 00351 void write_unformatted (std::ostream &out, const bool newline = true) const; 00352 00353 // protected: 00354 00358 T _coords[LIBMESH_DIM*LIBMESH_DIM]; 00359 }; 00360 00361 00362 template <typename T> 00363 class TypeTensorColumn 00364 { 00365 public: 00366 00367 TypeTensorColumn(TypeTensor<T> &tensor, 00368 unsigned int j) : 00369 _tensor(&tensor), _j(j) {} 00370 00375 T & operator () (const unsigned int i) 00376 { return (*_tensor)(i,_j); } 00377 00378 T & slice (const unsigned int i) 00379 { return (*_tensor)(i,_j); } 00380 00384 TypeTensorColumn<T>& operator = (const TypeVector<T>& rhs) 00385 { 00386 for (unsigned int i=0; i != LIBMESH_DIM; ++i) 00387 (*this)(i) = rhs(i); 00388 return *this; 00389 } 00390 00391 private: 00392 TypeTensor<T> *_tensor; 00393 const unsigned int _j; 00394 }; 00395 00396 00397 template <typename T> 00398 class ConstTypeTensorColumn 00399 { 00400 public: 00401 00402 ConstTypeTensorColumn(const TypeTensor<T> &tensor, 00403 unsigned int j) : 00404 _tensor(&tensor), _j(j) {} 00405 00409 const T & operator () (const unsigned int i) const 00410 { return (*_tensor)(i,_j); } 00411 00412 const T & slice (const unsigned int i) const 00413 { return (*_tensor)(i,_j); } 00414 00415 private: 00416 const TypeTensor<T> *_tensor; 00417 const unsigned int _j; 00418 }; 00419 00420 //------------------------------------------------------ 00421 // Inline functions 00422 template <typename T> 00423 inline 00424 TypeTensor<T>::TypeTensor () 00425 { 00426 _coords[0] = 0; 00427 00428 #if LIBMESH_DIM > 1 00429 _coords[1] = 0; 00430 _coords[2] = 0; 00431 _coords[3] = 0; 00432 #endif 00433 00434 #if LIBMESH_DIM > 2 00435 _coords[4] = 0; 00436 _coords[5] = 0; 00437 _coords[6] = 0; 00438 _coords[7] = 0; 00439 _coords[8] = 0; 00440 #endif 00441 } 00442 00443 00444 00445 template <typename T> 00446 inline 00447 TypeTensor<T>::TypeTensor 00448 (const T xx, 00449 const T xy, 00450 const T xz, 00451 const T yx, 00452 const T yy, 00453 const T yz, 00454 const T zx, 00455 const T zy, 00456 const T zz) 00457 { 00458 _coords[0] = xx; 00459 00460 #if LIBMESH_DIM == 2 00461 _coords[1] = xy; 00462 _coords[2] = yx; 00463 _coords[3] = yy; 00464 #endif 00465 00466 #if LIBMESH_DIM == 3 00467 _coords[1] = xy; 00468 _coords[2] = xz; 00469 _coords[3] = yx; 00470 _coords[4] = yy; 00471 _coords[5] = yz; 00472 _coords[6] = zx; 00473 _coords[7] = zy; 00474 _coords[8] = zz; 00475 #endif 00476 } 00477 00478 00479 template <typename T> 00480 template <typename Scalar> 00481 inline 00482 TypeTensor<T>::TypeTensor 00483 (const Scalar xx, 00484 const Scalar xy, 00485 const Scalar xz, 00486 const Scalar yx, 00487 const Scalar yy, 00488 const Scalar yz, 00489 const Scalar zx, 00490 const Scalar zy, 00491 typename 00492 boostcopy::enable_if_c<ScalarTraits<Scalar>::value, 00493 const Scalar>::type zz) 00494 { 00495 _coords[0] = xx; 00496 00497 #if LIBMESH_DIM == 2 00498 _coords[1] = xy; 00499 _coords[2] = yx; 00500 _coords[3] = yy; 00501 #endif 00502 00503 #if LIBMESH_DIM == 3 00504 _coords[1] = xy; 00505 _coords[2] = xz; 00506 _coords[3] = yx; 00507 _coords[4] = yy; 00508 _coords[5] = yz; 00509 _coords[6] = zx; 00510 _coords[7] = zy; 00511 _coords[8] = zz; 00512 #endif 00513 } 00514 00515 00516 00517 template <typename T> 00518 template<typename T2> 00519 inline 00520 TypeTensor<T>::TypeTensor (const TypeTensor <T2> &p) 00521 { 00522 // copy the nodes from vector p to me 00523 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00524 _coords[i] = p._coords[i]; 00525 } 00526 00527 00528 template <typename T> 00529 template <typename T2> 00530 TypeTensor<T>::TypeTensor(const TypeVector<T2>& vx) 00531 { 00532 libmesh_assert_equal_to (LIBMESH_DIM, 1); 00533 _coords[0] = vx(0); 00534 } 00535 00536 template <typename T> 00537 template <typename T2> 00538 TypeTensor<T>::TypeTensor(const TypeVector<T2>& vx, const TypeVector<T2> &vy) 00539 { 00540 libmesh_assert_equal_to (LIBMESH_DIM, 2); 00541 _coords[0] = vx(0); 00542 _coords[1] = vx(1); 00543 _coords[2] = vy(0); 00544 _coords[3] = vy(1); 00545 } 00546 00547 template <typename T> 00548 template <typename T2> 00549 TypeTensor<T>::TypeTensor(const TypeVector<T2>& vx, const TypeVector<T2> &vy, const TypeVector<T2> &vz) 00550 { 00551 libmesh_assert_equal_to (LIBMESH_DIM, 3); 00552 _coords[0] = vx(0); 00553 _coords[1] = vx(1); 00554 _coords[2] = vx(2); 00555 _coords[3] = vy(0); 00556 _coords[4] = vy(1); 00557 _coords[5] = vy(2); 00558 _coords[6] = vz(0); 00559 _coords[7] = vz(1); 00560 _coords[8] = vz(2); 00561 } 00562 00563 00564 00565 00566 template <typename T> 00567 inline 00568 TypeTensor<T>::~TypeTensor () 00569 { 00570 } 00571 00572 00573 00574 template <typename T> 00575 template<typename T2> 00576 inline 00577 void TypeTensor<T>::assign (const TypeTensor<T2> &p) 00578 { 00579 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00580 _coords[i] = p._coords[i]; 00581 } 00582 00583 00584 00585 template <typename T> 00586 inline 00587 const T & TypeTensor<T>::operator () (const unsigned int i, 00588 const unsigned int j) const 00589 { 00590 libmesh_assert_less (i, 3); 00591 libmesh_assert_less (j, 3); 00592 00593 #if LIBMESH_DIM < 3 00594 const static T my_zero = 0; 00595 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM) 00596 return my_zero; 00597 #endif 00598 00599 return _coords[i*LIBMESH_DIM+j]; 00600 } 00601 00602 00603 00604 template <typename T> 00605 inline 00606 T & TypeTensor<T>::operator () (const unsigned int i, 00607 const unsigned int j) 00608 { 00609 #if LIBMESH_DIM < 3 00610 00611 if (i >= LIBMESH_DIM || j >= LIBMESH_DIM) 00612 libmesh_error_msg("ERROR: You are assigning to a tensor component that is out of range for the compiled LIBMESH_DIM!"); 00613 00614 #endif 00615 00616 libmesh_assert_less (i, LIBMESH_DIM); 00617 libmesh_assert_less (j, LIBMESH_DIM); 00618 00619 return _coords[i*LIBMESH_DIM+j]; 00620 } 00621 00622 00623 template <typename T> 00624 inline 00625 ConstTypeTensorColumn<T> 00626 TypeTensor<T>::slice (const unsigned int i) const 00627 { 00628 libmesh_assert_less (i, LIBMESH_DIM); 00629 return ConstTypeTensorColumn<T>(*this, i); 00630 } 00631 00632 00633 template <typename T> 00634 inline 00635 TypeTensorColumn<T> 00636 TypeTensor<T>::slice (const unsigned int i) 00637 { 00638 libmesh_assert_less (i, LIBMESH_DIM); 00639 return TypeTensorColumn<T>(*this, i); 00640 } 00641 00642 00643 template <typename T> 00644 inline 00645 TypeVector<T> 00646 TypeTensor<T>::row(const unsigned int r) 00647 { 00648 TypeVector<T> return_vector; 00649 00650 for(unsigned int j=0; j<LIBMESH_DIM; j++) 00651 return_vector._coords[j] = _coords[r*LIBMESH_DIM + j]; 00652 00653 return return_vector; 00654 } 00655 00656 template <typename T> 00657 template<typename T2> 00658 inline 00659 TypeTensor<typename CompareTypes<T, T2>::supertype> 00660 TypeTensor<T>::operator + (const TypeTensor<T2> &p) const 00661 { 00662 00663 #if LIBMESH_DIM == 1 00664 return TypeTensor(_coords[0] + p._coords[0]); 00665 #endif 00666 00667 #if LIBMESH_DIM == 2 00668 return TypeTensor(_coords[0] + p._coords[0], 00669 _coords[1] + p._coords[1], 00670 0., 00671 _coords[2] + p._coords[2], 00672 _coords[3] + p._coords[3]); 00673 #endif 00674 00675 #if LIBMESH_DIM == 3 00676 return TypeTensor(_coords[0] + p._coords[0], 00677 _coords[1] + p._coords[1], 00678 _coords[2] + p._coords[2], 00679 _coords[3] + p._coords[3], 00680 _coords[4] + p._coords[4], 00681 _coords[5] + p._coords[5], 00682 _coords[6] + p._coords[6], 00683 _coords[7] + p._coords[7], 00684 _coords[8] + p._coords[8]); 00685 #endif 00686 00687 } 00688 00689 00690 00691 template <typename T> 00692 template<typename T2> 00693 inline 00694 const TypeTensor<T> & TypeTensor<T>::operator += (const TypeTensor<T2> &p) 00695 { 00696 this->add (p); 00697 00698 return *this; 00699 } 00700 00701 00702 00703 template <typename T> 00704 template<typename T2> 00705 inline 00706 void TypeTensor<T>::add (const TypeTensor<T2> &p) 00707 { 00708 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00709 _coords[i] += p._coords[i]; 00710 } 00711 00712 00713 00714 template <typename T> 00715 template <typename T2> 00716 inline 00717 void TypeTensor<T>::add_scaled (const TypeTensor<T2> &p, const T factor) 00718 { 00719 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00720 _coords[i] += factor*p._coords[i]; 00721 00722 } 00723 00724 00725 00726 template <typename T> 00727 template<typename T2> 00728 inline 00729 TypeTensor<typename CompareTypes<T, T2>::supertype> 00730 TypeTensor<T>::operator - (const TypeTensor<T2> &p) const 00731 { 00732 00733 #if LIBMESH_DIM == 1 00734 return TypeTensor(_coords[0] - p._coords[0]); 00735 #endif 00736 00737 #if LIBMESH_DIM == 2 00738 return TypeTensor(_coords[0] - p._coords[0], 00739 _coords[1] - p._coords[1], 00740 0., 00741 _coords[2] - p._coords[2], 00742 _coords[3] - p._coords[3]); 00743 #endif 00744 00745 #if LIBMESH_DIM == 3 00746 return TypeTensor(_coords[0] - p._coords[0], 00747 _coords[1] - p._coords[1], 00748 _coords[2] - p._coords[2], 00749 _coords[3] - p._coords[3], 00750 _coords[4] - p._coords[4], 00751 _coords[5] - p._coords[5], 00752 _coords[6] - p._coords[6], 00753 _coords[7] - p._coords[7], 00754 _coords[8] - p._coords[8]); 00755 #endif 00756 00757 } 00758 00759 00760 00761 template <typename T> 00762 template <typename T2> 00763 inline 00764 const TypeTensor<T> & TypeTensor<T>::operator -= (const TypeTensor<T2> &p) 00765 { 00766 this->subtract (p); 00767 00768 return *this; 00769 } 00770 00771 00772 00773 template <typename T> 00774 template <typename T2> 00775 inline 00776 void TypeTensor<T>::subtract (const TypeTensor<T2>& p) 00777 { 00778 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00779 _coords[i] -= p._coords[i]; 00780 } 00781 00782 00783 00784 template <typename T> 00785 template <typename T2> 00786 inline 00787 void TypeTensor<T>::subtract_scaled (const TypeTensor<T2> &p, const T factor) 00788 { 00789 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00790 _coords[i] -= factor*p._coords[i]; 00791 } 00792 00793 00794 00795 template <typename T> 00796 inline 00797 TypeTensor<T> TypeTensor<T>::operator - () const 00798 { 00799 00800 #if LIBMESH_DIM == 1 00801 return TypeTensor(-_coords[0]); 00802 #endif 00803 00804 #if LIBMESH_DIM == 2 00805 return TypeTensor(-_coords[0], 00806 -_coords[1], 00807 -_coords[2], 00808 -_coords[3]); 00809 #endif 00810 00811 #if LIBMESH_DIM == 3 00812 return TypeTensor(-_coords[0], 00813 -_coords[1], 00814 -_coords[2], 00815 -_coords[3], 00816 -_coords[4], 00817 -_coords[5], 00818 -_coords[6], 00819 -_coords[7], 00820 -_coords[8]); 00821 #endif 00822 00823 } 00824 00825 00826 00827 template <typename T> 00828 template <typename Scalar> 00829 inline 00830 typename boostcopy::enable_if_c< 00831 ScalarTraits<Scalar>::value, 00832 TypeTensor<typename CompareTypes<T, Scalar>::supertype> >::type 00833 TypeTensor<T>::operator * (const Scalar factor) const 00834 { 00835 typedef typename CompareTypes<T, Scalar>::supertype TS; 00836 00837 00838 #if LIBMESH_DIM == 1 00839 return TypeTensor<TS>(_coords[0]*factor); 00840 #endif 00841 00842 #if LIBMESH_DIM == 2 00843 return TypeTensor<TS>(_coords[0]*factor, 00844 _coords[1]*factor, 00845 _coords[2]*factor, 00846 _coords[3]*factor); 00847 #endif 00848 00849 #if LIBMESH_DIM == 3 00850 return TypeTensor<TS>(_coords[0]*factor, 00851 _coords[1]*factor, 00852 _coords[2]*factor, 00853 _coords[3]*factor, 00854 _coords[4]*factor, 00855 _coords[5]*factor, 00856 _coords[6]*factor, 00857 _coords[7]*factor, 00858 _coords[8]*factor); 00859 #endif 00860 } 00861 00862 00863 00864 template <typename T, typename Scalar> 00865 inline 00866 typename boostcopy::enable_if_c< 00867 ScalarTraits<Scalar>::value, 00868 TypeTensor<typename CompareTypes<T, Scalar>::supertype> >::type 00869 operator * (const Scalar factor, 00870 const TypeTensor<T> &t) 00871 { 00872 return t * factor; 00873 } 00874 00875 00876 00877 template <typename T> 00878 template <typename Scalar> 00879 inline 00880 const TypeTensor<T> & TypeTensor<T>::operator *= (const Scalar factor) 00881 { 00882 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 00883 _coords[i] *= factor; 00884 00885 return *this; 00886 } 00887 00888 00889 00890 00891 template <typename T> 00892 template <typename Scalar> 00893 inline 00894 typename boostcopy::enable_if_c< 00895 ScalarTraits<Scalar>::value, 00896 TypeTensor<typename CompareTypes<T, Scalar>::supertype> >::type 00897 TypeTensor<T>::operator / (const Scalar factor) const 00898 { 00899 libmesh_assert_not_equal_to (factor, static_cast<T>(0.)); 00900 00901 typedef typename CompareTypes<T, Scalar>::supertype TS; 00902 00903 #if LIBMESH_DIM == 1 00904 return TypeTensor<TS>(_coords[0]/factor); 00905 #endif 00906 00907 #if LIBMESH_DIM == 2 00908 return TypeTensor<TS>(_coords[0]/factor, 00909 _coords[1]/factor, 00910 _coords[2]/factor, 00911 _coords[3]/factor); 00912 #endif 00913 00914 #if LIBMESH_DIM == 3 00915 return TypeTensor<TS>(_coords[0]/factor, 00916 _coords[1]/factor, 00917 _coords[2]/factor, 00918 _coords[3]/factor, 00919 _coords[4]/factor, 00920 _coords[5]/factor, 00921 _coords[6]/factor, 00922 _coords[7]/factor, 00923 _coords[8]/factor); 00924 #endif 00925 00926 } 00927 00928 00929 00930 template <typename T> 00931 inline 00932 TypeTensor<T> TypeTensor<T>::transpose() const 00933 { 00934 #if LIBMESH_DIM == 1 00935 return TypeTensor(_coords[0]); 00936 #endif 00937 00938 #if LIBMESH_DIM == 2 00939 return TypeTensor(_coords[0], 00940 _coords[2], 00941 _coords[1], 00942 _coords[3]); 00943 #endif 00944 00945 #if LIBMESH_DIM == 3 00946 return TypeTensor(_coords[0], 00947 _coords[3], 00948 _coords[6], 00949 _coords[1], 00950 _coords[4], 00951 _coords[7], 00952 _coords[2], 00953 _coords[5], 00954 _coords[8]); 00955 #endif 00956 } 00957 00958 00959 00960 template <typename T> 00961 inline 00962 TypeTensor<T> TypeTensor<T>::inverse() const 00963 { 00964 #if LIBMESH_DIM == 1 00965 libmesh_assert_not_equal_to(_coords[0], static_cast<T>(0.)); 00966 return TypeTensor(1. / _coords[0]); 00967 #endif 00968 00969 #if LIBMESH_DIM == 2 00970 // Get convenient reference to this. 00971 const TypeTensor<T> & A = *this; 00972 00973 // Use temporary variables, avoid multiple accesses 00974 T a = A(0,0), b = A(0,1), 00975 c = A(1,0), d = A(1,1); 00976 00977 // Make sure det = ad - bc is not zero 00978 T my_det = a*d - b*c; 00979 libmesh_assert_not_equal_to(my_det, static_cast<T>(0.)); 00980 00981 return TypeTensor(d/my_det, -b/my_det, -c/my_det, a/my_det); 00982 #endif 00983 00984 #if LIBMESH_DIM == 3 00985 // Get convenient reference to this. 00986 const TypeTensor<T> & A = *this; 00987 00988 T a11 = A(0,0), a12 = A(0,1), a13 = A(0,2), 00989 a21 = A(1,0), a22 = A(1,1), a23 = A(1,2), 00990 a31 = A(2,0), a32 = A(2,1), a33 = A(2,2); 00991 00992 T my_det = a11*(a33*a22-a32*a23) - a21*(a33*a12-a32*a13) + a31*(a23*a12-a22*a13); 00993 libmesh_assert_not_equal_to(my_det, static_cast<T>(0.)); 00994 00995 // Inline comment characters are for lining up columns. 00996 return TypeTensor( (a33*a22-a32*a23)/my_det, -(a33*a12-a32*a13)/my_det, (a23*a12-a22*a13)/my_det, 00997 -(a33*a21-a31*a23)/my_det, (a33*a11-a31*a13)/my_det, -(a23*a11-a21*a13)/my_det, 00998 (a32*a21-a31*a22)/my_det, -(a32*a11-a31*a12)/my_det, (a22*a11-a21*a12)/my_det); 00999 #endif 01000 } 01001 01002 01003 01004 template <typename T> 01005 inline 01006 const TypeTensor<T> & TypeTensor<T>::operator /= (const T factor) 01007 { 01008 libmesh_assert_not_equal_to (factor, static_cast<T>(0.)); 01009 01010 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01011 _coords[i] /= factor; 01012 01013 return *this; 01014 } 01015 01016 01017 01018 01019 template <typename T> 01020 template <typename T2> 01021 inline 01022 TypeVector<typename CompareTypes<T,T2>::supertype> 01023 TypeTensor<T>::operator * (const TypeVector<T2> &p) const 01024 { 01025 TypeVector<typename CompareTypes<T,T2>::supertype> returnval; 01026 for (unsigned int i=0; i<LIBMESH_DIM; i++) 01027 for (unsigned int j=0; j<LIBMESH_DIM; j++) 01028 returnval(i) += (*this)(i,j)*p(j); 01029 01030 return returnval; 01031 } 01032 01033 01034 01035 template <typename T> 01036 template <typename T2> 01037 inline 01038 TypeTensor<T> TypeTensor<T>::operator * (const TypeTensor<T2> &p) const 01039 { 01040 TypeTensor<T> returnval; 01041 for (unsigned int i=0; i<LIBMESH_DIM; i++) 01042 for (unsigned int j=0; j<LIBMESH_DIM; j++) 01043 for (unsigned int k=0; k<LIBMESH_DIM; k++) 01044 returnval(i,j) += (*this)(i,k)*p(k,j); 01045 01046 return returnval; 01047 } 01048 01049 01050 01055 template <typename T> 01056 template <typename T2> 01057 inline 01058 typename CompareTypes<T,T2>::supertype 01059 TypeTensor<T>::contract (const TypeTensor<T2> &t) const 01060 { 01061 typename CompareTypes<T,T2>::supertype sum = 0.; 01062 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01063 sum += _coords[i]*t._coords[i]; 01064 return sum; 01065 } 01066 01067 01068 01069 template <typename T> 01070 inline 01071 Real TypeTensor<T>::size() const 01072 { 01073 return std::sqrt(this->size_sq()); 01074 } 01075 01076 01077 01078 template <typename T> 01079 inline 01080 T TypeTensor<T>::det() const 01081 { 01082 #if LIBMESH_DIM == 1 01083 return _coords[0]; 01084 #endif 01085 01086 #if LIBMESH_DIM == 2 01087 return (_coords[0] * _coords[3] 01088 - _coords[1] * _coords[2]); 01089 #endif 01090 01091 #if LIBMESH_DIM == 3 01092 return (_coords[0] * _coords[4] * _coords[8] 01093 + _coords[1] * _coords[5] * _coords[6] 01094 + _coords[2] * _coords[3] * _coords[7] 01095 - _coords[0] * _coords[5] * _coords[7] 01096 - _coords[1] * _coords[3] * _coords[8] 01097 - _coords[2] * _coords[4] * _coords[6]); 01098 #endif 01099 } 01100 01101 template <typename T> 01102 inline 01103 T TypeTensor<T>::tr() const 01104 { 01105 #if LIBMESH_DIM == 1 01106 return _coords[0]; 01107 #endif 01108 01109 #if LIBMESH_DIM == 2 01110 return _coords[0] + _coords[3]; 01111 #endif 01112 01113 #if LIBMESH_DIM == 3 01114 return _coords[0] + _coords[4] + _coords[8]; 01115 #endif 01116 } 01117 01118 template <typename T> 01119 inline 01120 void TypeTensor<T>::zero() 01121 { 01122 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01123 _coords[i] = 0.; 01124 } 01125 01126 01127 01128 template <typename T> 01129 inline 01130 Real TypeTensor<T>::size_sq () const 01131 { 01132 Real sum = 0.; 01133 for (unsigned int i=0; i<LIBMESH_DIM*LIBMESH_DIM; i++) 01134 sum += TensorTools::norm_sq(_coords[i]); 01135 return sum; 01136 } 01137 01138 01139 01140 template <typename T> 01141 inline 01142 bool TypeTensor<T>::operator == (const TypeTensor<T>& rhs) const 01143 { 01144 #if LIBMESH_DIM == 1 01145 return (std::abs(_coords[0] - rhs._coords[0]) 01146 < TOLERANCE); 01147 #endif 01148 01149 #if LIBMESH_DIM == 2 01150 return ((std::abs(_coords[0] - rhs._coords[0]) + 01151 std::abs(_coords[1] - rhs._coords[1]) + 01152 std::abs(_coords[2] - rhs._coords[2]) + 01153 std::abs(_coords[3] - rhs._coords[3])) 01154 < 4.*TOLERANCE); 01155 #endif 01156 01157 #if LIBMESH_DIM == 3 01158 return ((std::abs(_coords[0] - rhs._coords[0]) + 01159 std::abs(_coords[1] - rhs._coords[1]) + 01160 std::abs(_coords[2] - rhs._coords[2]) + 01161 std::abs(_coords[3] - rhs._coords[3]) + 01162 std::abs(_coords[4] - rhs._coords[4]) + 01163 std::abs(_coords[5] - rhs._coords[5]) + 01164 std::abs(_coords[6] - rhs._coords[6]) + 01165 std::abs(_coords[7] - rhs._coords[7]) + 01166 std::abs(_coords[8] - rhs._coords[8])) 01167 < 9.*TOLERANCE); 01168 #endif 01169 01170 } 01171 01172 } // namespace libMesh 01173 01174 #endif // LIBMESH_TYPE_TENSOR_H