void NeoHookeanMaterialNonLinear<Mesh>::computeStiffness ( const vector_Type&       sol,
                                                           Real                     /*factor*/,
                                                           const dataPtr_Type&      dataMaterial,
                                                           const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
                                                           const displayerPtr_Type& displayer )
{
    this->M_stiff.reset (new vector_Type (*this->M_localMap) );

    displayer->leaderPrint (" \n******************************************************************\n  ");
    displayer->leaderPrint (" Non-Linear S-  Computing the Neo-Hookean nonlinear stiffness vector"     );
    displayer->leaderPrint (" \n******************************************************************\n  ");

    UInt totalDof   = this->M_FESpace->dof().numTotalDof();
    UInt dim        = this->M_FESpace->dim();

    VectorElemental dk_loc ( this->M_FESpace->fe().nbFEDof(), nDimensions );
    //vector_Type disp(sol);

    vector_Type dRep (sol, Repeated);

    mapIterator_Type it;

    for ( it = (*mapsMarkerVolumes).begin(); it != (*mapsMarkerVolumes).end(); it++ )
    {

        //Given the marker pointed by the iterator, let's extract the material parameters
        UInt marker = it->first;

        Real mu     = dataMaterial->mu (marker);
        Real bulk   = dataMaterial->bulk (marker);

        for ( UInt j (0); j < it->second.size(); j++ )
        {
            this->M_FESpace->fe().updateFirstDerivQuadPt ( * (it->second[j]) );

            UInt eleID = this->M_FESpace->fe().currentLocalId();

            for ( UInt iNode = 0 ; iNode < ( UInt ) this->M_FESpace->fe().nbFEDof() ; iNode++ )
            {
                UInt  iloc = this->M_FESpace->fe().patternFirst ( iNode );

                for ( UInt iComp = 0; iComp < nDimensions; ++iComp )
                {
                    UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + iComp * dim + this->M_offset;
                    dk_loc[ iloc + iComp * this->M_FESpace->fe().nbFEDof() ] = dRep[ig];
                }
            }

            this->M_elvecK->zero();

            computeKinematicsVariables ( dk_loc );

            //! Stiffness for non-linear terms of the Neo-Hookean model
            /*!
              The results of the integrals are stored at each step into elvecK, until to build K matrix of the bilinear form
            */
            //! Volumetric part
            /*!
              Source term Pvol: int { bulk /2* (J1^2 - J1  + log(J1) ) * 1/J1 * (CofF1 : \nabla v) }
            */
            AssemblyElementalStructure::source_Pvol ( 0.5 * bulk, (*M_CofFk), (*M_Jack),
                                                      *this->M_elvecK,  this->M_FESpace->fe() );

            //! Isochoric part
            /*!
              Source term P1iso_NH
            */
            AssemblyElementalStructure::source_P1iso_NH ( mu, (*M_CofFk) , (*M_Fk),  (*M_Jack),  (*M_trCisok) ,
                                                          *this->M_elvecK,  this->M_FESpace->fe() );

            for ( UInt ic = 0; ic < nDimensions; ++ic )
            {
                /*!
                  M_elvecK is assemble into *vec_stiff vector that is recall
                  from updateSystem(matrix_ptrtype& mat_stiff, vector_ptr_type& vec_stiff)
                */
                assembleVector ( *this->M_stiff,
                                 *this->M_elvecK,
                                 this->M_FESpace->fe(),
                                 this->M_FESpace->dof(),
                                 ic, this->M_offset +  ic * totalDof );
            }
        }
    }

    this->M_stiff->globalAssemble();
}
void NeoHookeanMaterialNonLinear<Mesh>::updateNonLinearJacobianTerms ( matrixPtr_Type&       jacobian,
                                                                       const vector_Type&    disp,
                                                                       const dataPtr_Type&   dataMaterial,
                                                                       const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
                                                                       const displayerPtr_Type&  displayer )
{
    displayer->leaderPrint ("   Non-Linear S-  updating non linear terms in the Jacobian Matrix (Neo-Hookean)");

    UInt totalDof = this->M_FESpace->dof().numTotalDof();
    VectorElemental dk_loc (this->M_FESpace->fe().nbFEDof(), nDimensions);

    vector_Type dRep (disp, Repeated);

    //! Number of displacement components
    UInt nc = nDimensions;

    //! Nonlinear part of jacobian
    //! loop on volumes: assembling source term

    mapIterator_Type it;

    for ( it = (*mapsMarkerVolumes).begin(); it != (*mapsMarkerVolumes).end(); it++ )
    {

        //Given the marker pointed by the iterator, let's extract the material parameters
        UInt marker = it->first;

        Real mu     = dataMaterial->mu (marker);
        Real bulk   = dataMaterial->bulk (marker);

        for ( UInt j (0); j < it->second.size(); j++ )
        {
            this->M_FESpace->fe().updateFirstDerivQuadPt ( * (it->second[j]) );

            UInt eleID = this->M_FESpace->fe().currentLocalId();

            for ( UInt iNode = 0; iNode < ( UInt ) this->M_FESpace->fe().nbFEDof(); iNode++ )
            {
                UInt  iloc = this->M_FESpace->fe().patternFirst ( iNode );

                for ( UInt iComp = 0; iComp < nDimensions; ++iComp )
                {
                    UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + iComp * this->M_FESpace->dim() + this->M_offset;
                    dk_loc[iloc + iComp * this->M_FESpace->fe().nbFEDof()] = dRep[ig];
                }
            }

            this->M_elmatK->zero();

            //! Computes F, Cof(F), J = det(F), Tr(C)
            computeKinematicsVariables ( dk_loc );

            //! Stiffness for non-linear terms of the Neo-Hookean model
            /*!
              The results of the integrals are stored at each step into elmatK, until to build K matrix of the bilinear form
            */

            //! VOLUMETRIC PART
            //! 1. Stiffness matrix: int { 1/2 * bulk * ( 2 - 1/J + 1/J^2 ) * ( CofF : \nabla \delta ) (CofF : \nabla v) }
            AssemblyElementalStructure::stiff_Jac_Pvol_1term ( 0.5 * bulk, (*M_CofFk), (*M_Jack),
                                                               *this->M_elmatK, this->M_FESpace->fe() );

            //! 2. Stiffness matrix: int { 1/2 * bulk * ( 1/J- 1 - log(J)/J^2 ) * ( CofF [\nabla \delta]^t CofF ) : \nabla v }
            AssemblyElementalStructure::stiff_Jac_Pvol_2term ( 0.5 * bulk, (*M_CofFk), (*M_Jack),
                                                               *this->M_elmatK, this->M_FESpace->fe() );

            //! ISOCHORIC PART
            //! 1. Stiffness matrix : int { -2/3 * mu * J^(-5/3) *( CofF : \nabla \delta ) ( F : \nabla \v ) }
            AssemblyElementalStructure::stiff_Jac_P1iso_NH_1term ( (-2.0 / 3.0) * mu, (*M_CofFk), (*M_Fk), (*M_Jack),
                                                                   *this->M_elmatK, this->M_FESpace->fe() );

            //! 2. Stiffness matrix : int { 2/9 * mu * ( Ic_iso / J^2 )( CofF : \nabla \delta ) ( CofF : \nabla \v ) }
            AssemblyElementalStructure::stiff_Jac_P1iso_NH_2term ( (2.0 / 9.0) * mu, (*M_CofFk), (*M_Jack), (*M_trCisok),
                                                                   *this->M_elmatK, this->M_FESpace->fe() );

            //! 3. Stiffness matrix : int { mu * J^(-2/3) (\nabla \delta : \nabla \v)}
            AssemblyElementalStructure::stiff_Jac_P1iso_NH_3term ( mu, (*M_Jack), *this->M_elmatK, this->M_FESpace->fe() );

            //! 4. Stiffness matrix : int { -2/3 * mu * J^(-5/3) ( F : \nabla \delta ) ( CofF : \nabla \v ) }
            AssemblyElementalStructure::stiff_Jac_P1iso_NH_4term ( (-2.0 / 3.0) * mu, (*M_CofFk), (*M_Fk), (*M_Jack),
                                                                   *this->M_elmatK, this->M_FESpace->fe() );

            //! 5. Stiffness matrix : int { 1/3 * mu * J^(-2) * Ic_iso * (CofF [\nabla \delta]^t CofF ) : \nabla \v }
            AssemblyElementalStructure::stiff_Jac_P1iso_NH_5term ( (1.0 / 3.0) * mu, (*M_CofFk), (*M_Jack), (*M_trCisok),
                                                                   *this->M_elmatK, this->M_FESpace->fe() );

            //! assembling
            for ( UInt ic = 0; ic < nc; ++ic )
            {
                for ( UInt jc = 0; jc < nc; jc++ )
                {
                    assembleMatrix ( *jacobian,
                                     *this->M_elmatK,
                                     this->M_FESpace->fe(),
                                     this->M_FESpace->dof(),
                                     ic, jc,
                                     this->M_offset +  ic * totalDof, this->M_offset +  jc * totalDof );
                }
            }
        }
    }
}
void
WallTensionEstimatorCylindricalCoordinates<Mesh >::analyzeTensionsRecoveryEigenvaluesCylindrical ( void )
{

    LifeChrono chrono;

    this->M_displayer->leaderPrint (" \n*********************************\n  ");
    this->M_displayer->leaderPrint ("   Performing the analysis recovering the tensions..., ", this->M_dataMaterial->solidType() );
    this->M_displayer->leaderPrint (" \n*********************************\n  ");

    solutionVect_Type patchArea (* (this->M_displacement), Unique, Add);
    patchArea *= 0.0;

    super::constructPatchAreaVector ( patchArea );

    //Before assembling the reconstruction process is done
    solutionVect_Type patchAreaR (patchArea, Repeated);

    QuadratureRule fakeQuadratureRule;

    Real refElemArea (0); //area of reference element
    //compute the area of reference element
    for (UInt iq = 0; iq < this->M_FESpace->qr().nbQuadPt(); iq++)
    {
        refElemArea += this->M_FESpace->qr().weight (iq);
    }

    Real wQuad (refElemArea / this->M_FESpace->refFE().nbDof() );

    //Setting the quadrature Points = DOFs of the element and weight = 1
    std::vector<GeoVector> coords = this->M_FESpace->refFE().refCoor();
    std::vector<Real> weights (this->M_FESpace->fe().nbFEDof(), wQuad);
    fakeQuadratureRule.setDimensionShape ( shapeDimension (this->M_FESpace->refFE().shape() ), this->M_FESpace->refFE().shape() );
    fakeQuadratureRule.setPoints (coords, weights);

    //Set the new quadrature rule
    this->M_FESpace->setQuadRule (fakeQuadratureRule);

    UInt totalDof = this->M_FESpace->dof().numTotalDof();
    VectorElemental dk_loc (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );

    //Vectors for the deformation tensor
    std::vector<matrix_Type> vectorDeformationF (this->M_FESpace->fe().nbFEDof(), * (this->M_deformationF) );
    //Copying the displacement field into a vector with repeated map for parallel computations
    solutionVect_Type dRep (* (this->M_displacement), Repeated);

    VectorElemental elVecTens (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );

    chrono.start();

    //Loop on each volume
    for ( UInt i = 0; i < this->M_FESpace->mesh()->numVolumes(); ++i )
    {
        this->M_FESpace->fe().updateFirstDerivQuadPt ( this->M_FESpace->mesh()->volumeList ( i ) );
        elVecTens.zero();

        this->M_marker = this->M_FESpace->mesh()->volumeList ( i ).markerID();

        UInt eleID = this->M_FESpace->fe().currentLocalId();

        //Extracting the local displacement
        for ( UInt iNode = 0; iNode < ( UInt ) this->M_FESpace->fe().nbFEDof(); iNode++ )
        {
            UInt  iloc = this->M_FESpace->fe().patternFirst ( iNode );

            for ( UInt iComp = 0; iComp < this->M_FESpace->fieldDim(); ++iComp )
            {
                UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + iComp * this->M_FESpace->dim() + this->M_offset;
                dk_loc[iloc + iComp * this->M_FESpace->fe().nbFEDof()] = dRep[ig];
            }
        }

        //Compute the element tensor F
        AssemblyElementalStructure::computeLocalDeformationGradientWithoutIdentity ( dk_loc, vectorDeformationF, this->M_FESpace->fe() );

        //Compute the local vector of the principal stresses
        for ( UInt nDOF = 0; nDOF < ( UInt ) this->M_FESpace->fe().nbFEDof(); nDOF++ )
        {
            UInt  iloc = this->M_FESpace->fe().patternFirst ( nDOF );
            vector_Type localDisplacement (this->M_FESpace->fieldDim(), 0.0);

            for ( UInt coor = 0; coor < this->M_FESpace->fieldDim(); coor++ )
            {
                localDisplacement[coor] = iloc + coor * this->M_FESpace->fe().nbFEDof();
            }

            this->M_sigma->Scale (0.0);
            this->M_firstPiola->Scale (0.0);
            this->M_cofactorF->Scale (0.0);
            M_deformationCylindricalF->Scale (0.0);

            moveToCylindricalCoordinates (vectorDeformationF[nDOF], iloc, *M_deformationCylindricalF);

            //Compute the rightCauchyC tensor
            AssemblyElementalStructure::computeInvariantsRightCauchyGreenTensor (this->M_invariants, *M_deformationCylindricalF, * (this->M_cofactorF) );

            //Compute the first Piola-Kirchhoff tensor
            this->M_material->computeLocalFirstPiolaKirchhoffTensor (* (this->M_firstPiola), *M_deformationCylindricalF, * (this->M_cofactorF), this->M_invariants, this->M_marker);

            //Compute the Cauchy tensor
            AssemblyElementalStructure::computeCauchyStressTensor (* (this->M_sigma), * (this->M_firstPiola), this->M_invariants[3], *M_deformationCylindricalF);

            //Compute the eigenvalue
            AssemblyElementalStructure::computeEigenvalues (* (this->M_sigma), this->M_eigenvaluesR, this->M_eigenvaluesI);

            //The Cauchy tensor is symmetric and therefore, the eigenvalues are real
            //Check on the imaginary part of eigen values given by the Lapack method
            Real sum (0);
            for ( int i = 0; i < this->M_eigenvaluesI.size(); i++ )
            {
                sum += std::abs (this->M_eigenvaluesI[i]);
            }
            ASSERT_PRE ( sum < 1e-6 , "The eigenvalues of the Cauchy stress tensors have to be real!" );

            std::sort ( this->M_eigenvaluesR.begin(), this->M_eigenvaluesR.end() );

            //Assembling the local vector
            for ( int coor = 0; coor < this->M_eigenvaluesR.size(); coor++ )
            {
                elVecTens[iloc + coor * this->M_FESpace->fe().nbFEDof()] = this->M_eigenvaluesR[coor];
            }
        }

        super::reconstructElementaryVector ( elVecTens, patchAreaR, *this->M_FESpace );

        //Assembling the local into global vector
        for ( UInt ic = 0; ic < this->M_FESpace->fieldDim(); ++ic )
        {
            assembleVector (* (this->M_globalEigenvalues), elVecTens, this->M_FESpace->fe(), this->M_FESpace->dof(), ic, this->M_offset +  ic * totalDof );
        }
    }

    this->M_globalEigenvalues->globalAssemble();

    chrono.stop();
    this->M_displayer->leaderPrint ("Analysis done in: ", chrono.diff() );
}
void
WallTensionEstimatorCylindricalCoordinates<Mesh >::constructGlobalStressVector()
{

    //Creating the local stress tensors
    VectorElemental elVecSigmaX (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );
    VectorElemental elVecSigmaY (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );
    VectorElemental elVecSigmaZ (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );

    LifeChrono chrono;

    //Constructing the patch area vector for reconstruction purposes
    solutionVect_Type patchArea (* (this->M_displacement), Unique, Add);
    patchArea *= 0.0;

    super::constructPatchAreaVector ( patchArea );

    //Before assembling the reconstruction process is done
    solutionVect_Type patchAreaR (patchArea, Repeated);

    QuadratureRule fakeQuadratureRule;

    Real refElemArea (0); //area of reference element
    //compute the area of reference element
    for (UInt iq = 0; iq < this->M_FESpace->qr().nbQuadPt(); iq++)
    {
        refElemArea += this->M_FESpace->qr().weight (iq);
    }

    Real wQuad (refElemArea / this->M_FESpace->refFE().nbDof() );

    //Setting the quadrature Points = DOFs of the element and weight = 1
    std::vector<GeoVector> coords = this->M_FESpace->refFE().refCoor();
    std::vector<Real> weights (this->M_FESpace->fe().nbFEDof(), wQuad);
    fakeQuadratureRule.setDimensionShape ( shapeDimension (this->M_FESpace->refFE().shape() ), this->M_FESpace->refFE().shape() );
    fakeQuadratureRule.setPoints (coords, weights);

    //Set the new quadrature rule
    this->M_FESpace->setQuadRule (fakeQuadratureRule);

    this->M_displayer->leaderPrint (" \n*********************************\n  ");
    this->M_displayer->leaderPrint ("   Performing the analysis recovering the Cauchy stresses..., ", this->M_dataMaterial->solidType() );
    this->M_displayer->leaderPrint (" \n*********************************\n  ");

    UInt totalDof = this->M_FESpace->dof().numTotalDof();
    VectorElemental dk_loc (this->M_FESpace->fe().nbFEDof(), this->M_FESpace->fieldDim() );

    //Vectors for the deformation tensor
    std::vector<matrix_Type> vectorDeformationF (this->M_FESpace->fe().nbFEDof(), * (this->M_deformationF) );
    //Copying the displacement field into a vector with repeated map for parallel computations
    solutionVect_Type dRep (* (this->M_displacement), Repeated);

    chrono.start();

    //Loop on each volume
    for ( UInt i = 0; i < this->M_FESpace->mesh()->numVolumes(); ++i )
    {
        this->M_FESpace->fe().updateFirstDerivQuadPt ( this->M_FESpace->mesh()->volumeList ( i ) );

        elVecSigmaX.zero();
        elVecSigmaY.zero();
        elVecSigmaZ.zero();

        this->M_marker = this->M_FESpace->mesh()->volumeList ( i ).markerID();

        UInt eleID = this->M_FESpace->fe().currentLocalId();

        //Extracting the local displacement
        for ( UInt iNode = 0; iNode < ( UInt ) this->M_FESpace->fe().nbFEDof(); iNode++ )
        {
            UInt  iloc = this->M_FESpace->fe().patternFirst ( iNode );

            for ( UInt iComp = 0; iComp < this->M_FESpace->fieldDim(); ++iComp )
            {
                UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + iComp * this->M_FESpace->dim() + this->M_offset;
                dk_loc[iloc + iComp * this->M_FESpace->fe().nbFEDof()] = dRep[ig];
            }
        }

        //Compute the element tensor F
        AssemblyElementalStructure::computeLocalDeformationGradientWithoutIdentity ( dk_loc, vectorDeformationF, this->M_FESpace->fe() );

        //Compute the local vector of the principal stresses
        for ( UInt nDOF = 0; nDOF < ( UInt ) this->M_FESpace->fe().nbFEDof(); nDOF++ )
        {
            UInt  iloc = this->M_FESpace->fe().patternFirst ( nDOF );

            vector_Type localDisplacement (this->M_FESpace->fieldDim(), 0.0);

            for ( UInt coor = 0; coor < this->M_FESpace->fieldDim(); coor++ )
            {
                localDisplacement[coor] = iloc + coor * this->M_FESpace->fe().nbFEDof();
            }

            this->M_sigma->Scale (0.0);
            this->M_firstPiola->Scale (0.0);
            this->M_cofactorF->Scale (0.0);
            this->M_deformationCylindricalF->Scale (0.0);

            moveToCylindricalCoordinates (vectorDeformationF[nDOF], iloc, *M_deformationCylindricalF);

            //Compute the rightCauchyC tensor
            AssemblyElementalStructure::computeInvariantsRightCauchyGreenTensor (this->M_invariants, *M_deformationCylindricalF, * (this->M_cofactorF) );

            //Compute the first Piola-Kirchhoff tensor
            this->M_material->computeLocalFirstPiolaKirchhoffTensor (* (this->M_firstPiola), *M_deformationCylindricalF, * (this->M_cofactorF), this->M_invariants, this->M_marker);

            //Compute the Cauchy tensor
            AssemblyElementalStructure::computeCauchyStressTensor (* (this->M_sigma), * (this->M_firstPiola), this->M_invariants[3], *M_deformationCylindricalF);

            //Assembling the local vectors for local tensions Component X
            for ( int coor = 0; coor < this->M_FESpace->fieldDim(); coor++ )
            {
                (elVecSigmaX) [iloc + coor * this->M_FESpace->fe().nbFEDof()] = (* (this->M_sigma) ) (coor, 0);
            }

            //Assembling the local vectors for local tensions Component Y
            for ( int coor = 0; coor < this->M_FESpace->fieldDim(); coor++ )
            {
                (elVecSigmaY) [iloc + coor * this->M_FESpace->fe().nbFEDof()] = (* (this->M_sigma) ) (coor, 1);
            }

            //Assembling the local vectors for local tensions Component Z
            for ( int coor = 0; coor < this->M_FESpace->fieldDim(); coor++ )
            {
                (elVecSigmaZ) [iloc + coor * this->M_FESpace->fe().nbFEDof()] = (* (this->M_sigma) ) (coor, 2);
            }

        }

        super::reconstructElementaryVector ( elVecSigmaX, patchAreaR, *this->M_FESpace );
        super::reconstructElementaryVector ( elVecSigmaY, patchAreaR, *this->M_FESpace );
        super::reconstructElementaryVector ( elVecSigmaZ, patchAreaR, *this->M_FESpace );

        //Assembling the three elemental vector in the three global
        for ( UInt ic = 0; ic < this->M_FESpace->fieldDim(); ++ic )
        {
            assembleVector (*this->M_sigmaX, elVecSigmaX, this->M_FESpace->fe(), this->M_FESpace->dof(), ic, this->M_offset +  ic * totalDof );
            assembleVector (*this->M_sigmaY, elVecSigmaY, this->M_FESpace->fe(), this->M_FESpace->dof(), ic, this->M_offset +  ic * totalDof );
            assembleVector (*this->M_sigmaZ, elVecSigmaZ, this->M_FESpace->fe(), this->M_FESpace->dof(), ic, this->M_offset +  ic * totalDof );
        }
    }


    this->M_sigmaX->globalAssemble();
    this->M_sigmaY->globalAssemble();
    this->M_sigmaZ->globalAssemble();
}
void
WallTensionEstimatorCylindricalCoordinates<Mesh >::analyzeTensionsRecoveryDisplacementCylindrical ( void )
{

    LifeChrono chrono;

    this->M_displayer->leaderPrint (" \n*********************************\n  ");
    this->M_displayer->leaderPrint ("   Performing the analysis recovering the displacement..., ", this->M_dataMaterial->solidType() );
    this->M_displayer->leaderPrint (" \n*********************************\n  ");

    solutionVectPtr_Type grDisplX ( new solutionVect_Type (* (this->M_FESpace->mapPtr() ) ) );
    solutionVectPtr_Type grDisplY ( new solutionVect_Type (* (this->M_FESpace->mapPtr() ) ) );
    solutionVectPtr_Type grDisplZ ( new solutionVect_Type (* (this->M_FESpace->mapPtr() ) ) );

    //Compute the deformation gradient tensor F of the displ field
    super::computeDisplacementGradient ( grDisplX, grDisplY, grDisplZ);

    solutionVect_Type grXRep (*grDisplX, Repeated);
    solutionVect_Type grYRep (*grDisplY, Repeated);
    solutionVect_Type grZRep (*grDisplZ, Repeated);
    solutionVect_Type dRep (* (this->M_displacement), Repeated);

    //For each of the DOF, the Cauchy tensor is computed.
    //Therefore the tensor C,P, \sigma are computed for each DOF

    chrono.start();

    for ( UInt i (0); i < this->M_FESpace->mesh()->numVolumes(); ++i )
    {

        //Setup quantities
        this->M_FESpace->fe().updateFirstDerivQuadPt ( this->M_FESpace->mesh()->volumeList ( i ) );
        UInt eleID = this->M_FESpace->fe().currentLocalId();
        this->M_marker = this->M_FESpace->mesh()->volumeList ( i ).markerID();

        //store the local \grad(u)
        matrix_Type gradientDispl ( this->M_FESpace->fieldDim(), this->M_FESpace->fieldDim() );
        gradientDispl.Scale ( 0.0 );

        //Extracting the local displacement
        for ( UInt iNode = 0; iNode < ( UInt ) this->M_FESpace->fe().nbFEDof(); iNode++ )
        {
            UInt  iloc = this->M_FESpace->fe().patternFirst ( iNode );

            for ( UInt iComp = 0; iComp < this->M_FESpace->fieldDim(); ++iComp )
            {
                UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + iComp * this->M_FESpace->dim() + this->M_offset;

                gradientDispl (iComp, 0) = grXRep[ig];
                gradientDispl (iComp, 1) = grYRep[ig];
                gradientDispl (iComp, 2) = grZRep[ig];
            }

            //Reinitialization of matrices and arrays
            ( * (this->M_cofactorF) ).Scale (0.0);
            ( * (this->M_firstPiola) ).Scale (0.0);
            ( * (this->M_sigma) ).Scale (0.0);

            //Moving to cylindrical coordinates
            moveToCylindricalCoordinates ( gradientDispl, iloc, *M_deformationCylindricalF );

            //Compute the rightCauchyC tensor
            AssemblyElementalStructure::computeInvariantsRightCauchyGreenTensor (this->M_invariants, *M_deformationCylindricalF, * (this->M_cofactorF) );

            //Compute the first Piola-Kirchhoff tensor
            this->M_material->computeLocalFirstPiolaKirchhoffTensor (* (this->M_firstPiola), *M_deformationCylindricalF, * (this->M_cofactorF), this->M_invariants, this->M_marker);

            //Compute the Cauchy tensor
            AssemblyElementalStructure::computeCauchyStressTensor (* (this->M_sigma), * (this->M_firstPiola), this->M_invariants[3], *M_deformationCylindricalF);

            //Compute the eigenvalue
            AssemblyElementalStructure::computeEigenvalues (* (this->M_sigma), this->M_eigenvaluesR, this->M_eigenvaluesI);

            //The Cauchy tensor is symmetric and therefore, the eigenvalues are real
            //Check on the imaginary part of eigen values given by the Lapack method
            Real sum (0);
            for ( UInt j (0); j < this->M_eigenvaluesI.size(); ++j )
            {
                sum += std::abs ( this->M_eigenvaluesI[j] );
            }

            ASSERT ( sum < 1e-6 , "The eigenvalues of the Cauchy stress tensors have to be real!" );

            std::sort ( this->M_eigenvaluesR.begin(), this->M_eigenvaluesR.end() );

            std::cout << "Saving eigenvalues" << i << std::endl;

            //Save the eigenvalues in the global vector
            for ( UInt icoor = 0; icoor < this->M_FESpace->fieldDim(); ++icoor )
            {
                UInt ig = this->M_FESpace->dof().localToGlobalMap ( eleID, iloc ) + icoor * this->M_FESpace->dim() + this->M_offset;
                (* (this->M_globalEigenvalues) ) (ig) = this->M_eigenvaluesR[icoor];
                // Int LIDid = this->M_displacement->blockMap().LID(ig);
                // if( this->M_globalEigenvalues->blockMap().LID(ig) != -1  )
                // {
                //     Int GIDid = this->M_globalEigenvalues->blockMap().GID(LIDid);

                // }
            }
        }
    }

    chrono.stop();
    this->M_displayer->leaderPrint ("Analysis done in: ", chrono.diff() );

}