Skip to content
Snippets Groups Projects
tensor.C 9.51 KiB
Newer Older
  • Learn to ignore specific revisions
  • /*---------------------------------------------------------------------------*\
      =========                 |
      \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
       \\    /   O peration     |
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        \\  /    A nd           | www.openfoam.com
    
         \\/     M anipulation  |
    
    -------------------------------------------------------------------------------
    
    OpenFOAM bot's avatar
    OpenFOAM bot committed
        Copyright (C) 2011-2017 OpenFOAM Foundation
    
        Copyright (C) 2020 OpenCFD Ltd.
    
    -------------------------------------------------------------------------------
    License
        This file is part of OpenFOAM.
    
    
        OpenFOAM is free software: you can redistribute it and/or modify it
        under the terms of the GNU General Public License as published by
        the Free Software Foundation, either version 3 of the License, or
        (at your option) any later version.
    
    
        OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
        ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
        FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
        for more details.
    
        You should have received a copy of the GNU General Public License
    
        along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
    
    
    \*---------------------------------------------------------------------------*/
    
    #include "tensor.H"
    
    andy's avatar
    andy committed
    using namespace Foam::constant::mathematical;
    
    // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
    
    template<>
    const char* const Foam::tensor::vsType::typeName = "tensor";
    
    template<>
    const char* const Foam::tensor::vsType::componentNames[] =
    
        "xx", "xy", "xz",
        "yx", "yy", "yz",
        "zx", "zy", "zz"
    };
    
    template<>
    const Foam::tensor Foam::tensor::vsType::zero(tensor::uniform(0));
    
    template<>
    const Foam::tensor Foam::tensor::vsType::one(tensor::uniform(1));
    
    template<>
    const Foam::tensor Foam::tensor::vsType::max(tensor::uniform(VGREAT));
    
    template<>
    const Foam::tensor Foam::tensor::vsType::min(tensor::uniform(-VGREAT));
    
    template<>
    const Foam::tensor Foam::tensor::vsType::rootMax(tensor::uniform(ROOTVGREAT));
    
    template<>
    const Foam::tensor Foam::tensor::vsType::rootMin(tensor::uniform(-ROOTVGREAT));
    
    template<>
    const Foam::tensor Foam::tensor::I
    (
        1, 0, 0,
        0, 1, 0,
        0, 0, 1
    );
    
    
    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
    
    
    Foam::Vector<Foam::complex> Foam::eigenValues(const tensor& T)
    
        // Return diagonal if T is effectively diagonal tensor
        if
        (
            (
                sqr(T.xy()) + sqr(T.xz()) + sqr(T.yz())
              + sqr(T.yx()) + sqr(T.zx()) + sqr(T.zy())
            ) < ROOTSMALL
        )
        {
            return Vector<complex>
            (
                complex(T.xx()), complex(T.yy()), complex(T.zz())
            );
        }
    
    
        // Coefficients of the characteristic cubic polynomial (a = 1)
    
        const scalar b = - T.xx() - T.yy() - T.zz();
    
            T.xx()*T.yy() + T.xx()*T.zz() + T.yy()*T.zz()
          - T.xy()*T.yx() - T.yz()*T.zy() - T.zx()*T.xz();
    
          - T.xx()*T.yy()*T.zz()
          - T.xy()*T.yz()*T.zx() - T.xz()*T.zy()*T.yx()
          + T.xx()*T.yz()*T.zy() + T.yy()*T.zx()*T.xz() + T.zz()*T.xy()*T.yx();
    
        // Determine the roots of the characteristic cubic polynomial
        const Roots<3> roots(cubicEqn(1, b, c, d).roots());
    
        bool isComplex = false;
    
                    break;
                case roots::posInf:
                case roots::negInf:
                case roots::nan:
    
                    WarningInFunction
                        << "Eigenvalue computation fails for tensor: " << T
    
                        << "due to the not-a-number root = " << roots[i]
    
                case roots::real:
                    break;
    
            return
                Vector<complex>
                (
                    complex(roots[0], 0),
                    complex(roots[1], roots[2]),
                    complex(roots[1], -roots[2])
                );
    
        return
            Vector<complex>
            (
                complex(roots[0], 0),
                complex(roots[1], 0),
                complex(roots[2], 0)
            );
    
    Foam::Vector<Foam::complex> Foam::eigenVector
    
        const complex& eVal,
        const Vector<complex>& standardBasis1,
        const Vector<complex>& standardBasis2
    
        // Construct the characteristic equation system for this eigenvalue
        Tensor<complex> A(Zero);
        forAll(A, i)
        {
            A[i] = complex(T[i], 0);
        }
        A.xx() -= eVal;
        A.yy() -= eVal;
        A.zz() -= eVal;
    
        // Determinants of the 2x2 sub-matrices used to find the eigenvectors
    
        // Sub-determinants for a unique eigenvenvalue
        complex sd0 = A.yy()*A.zz() - A.yz()*A.zy();
        complex sd1 = A.zz()*A.xx() - A.zx()*A.xz();
        complex sd2 = A.xx()*A.yy() - A.xy()*A.yx();
        scalar magSd0 = mag(sd0);
        scalar magSd1 = mag(sd1);
        scalar magSd2 = mag(sd2);
    
    
        // Evaluate the eigenvector using the largest sub-determinant
    
        if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
    
            const Vector<complex> eVec
    
                (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
                (A.zy()*A.yx() - A.yy()*A.zx())/sd0
            );
    
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
            const Vector<complex> eVec
    
            (
                (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
    
                (A.zx()*A.xy() - A.xx()*A.zy())/sd1
            );
    
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
        }
        else if (magSd2 > SMALL)
        {
    
            const Vector<complex> eVec
    
            (
                (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
                (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
        // Sub-determinants for a repeated eigenvalue
    
        sd0 = A.yy()*standardBasis1.z() - A.yz()*standardBasis1.y();
        sd1 = A.zz()*standardBasis1.x() - A.zx()*standardBasis1.z();
        sd2 = A.xx()*standardBasis1.y() - A.xy()*standardBasis1.x();
    
        magSd0 = mag(sd0);
        magSd1 = mag(sd1);
        magSd2 = mag(sd2);
    
    
        // Evaluate the eigenvector using the largest sub-determinant
    
        if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
    
            const Vector<complex> eVec
    
                complex(1, 0),
                (A.yz()*standardBasis1.x() - standardBasis1.z()*A.yx())/sd0,
                (standardBasis1.y()*A.yx() - A.yy()*standardBasis1.x())/sd0
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
            const Vector<complex> eVec
    
                (standardBasis1.z()*A.zy() - A.zz()*standardBasis1.y())/sd1,
                complex(1, 0),
                (A.zx()*standardBasis1.y() - standardBasis1.x()*A.zy())/sd1
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
        }
        else if (magSd2 > SMALL)
        {
    
            const Vector<complex> eVec
    
                (A.xy()*standardBasis1.z() - standardBasis1.y()*A.xz())/sd2,
                (standardBasis1.x()*A.xz() - A.xx()*standardBasis1.z())/sd2,
                complex(1, 0)
    
            #ifdef FULLDEBUG
            if (mag(eVec) < SMALL)
            {
                FatalErrorInFunction
                    << "Eigenvector magnitude should be non-zero:"
                    << "mag(eigenvector) = " << mag(eVec)
                    << abort(FatalError);
            }
            #endif
    
            return eVec/mag(eVec);
    
        return standardBasis1^standardBasis2;
    
    Foam::Tensor<Foam::complex> Foam::eigenVectors
    
        const tensor& T,
        const Vector<complex>& eVals
    
        Vector<complex> Ux(complex(1, 0), Zero, Zero);
        Vector<complex> Uy(Zero, complex(1, 0), Zero);
        Vector<complex> Uz(Zero, Zero, complex(1, 0));
    
        Ux = eigenVector(T, eVals.x(), Uy, Uz);
        Uy = eigenVector(T, eVals.y(), Uz, Ux);
        Uz = eigenVector(T, eVals.z(), Ux, Uy);
    
        return Tensor<complex>(Ux, Uy, Uz);
    
    Foam::Tensor<Foam::complex> Foam::eigenVectors(const tensor& T)
    
        const Vector<complex> eVals(eigenValues(T));
    
        return eigenVectors(T, eVals);
    
    // ************************************************************************* //