diff --git a/src/cpp/tardigrade_hydra.cpp b/src/cpp/tardigrade_hydra.cpp index aa92271..70e2c82 100644 --- a/src/cpp/tardigrade_hydra.cpp +++ b/src/cpp/tardigrade_hydra.cpp @@ -1420,68 +1420,79 @@ namespace tardigradeHydra{ } - void hydraBase::solveNonLinearProblem( ){ + void hydraBase::solveNewtonUpdate( floatVector &deltaX_tr ){ /*! - * Solve the non-linear problem + * Solve the Newton update returning the trial value of the unknown vector + * + * \param &deltaX_tr: The trial change in the unknown vector */ - // Form the initial unknown vector - TARDIGRADE_ERROR_TOOLS_CATCH( initializeUnknownVector( ) ); - unsigned int rank; - floatVector deltaX( getNumUnknowns( ), 0 ); + Eigen::Map< Eigen::Vector< floatType, -1 > > dx_map( deltaX_tr.data( ), getNumUnknowns( ) ); - Eigen::Map< Eigen::Vector< floatType, -1 > > dx_map( deltaX.data( ), getUnknownVector( )->size( ) ); + tardigradeVectorTools::solverType< floatType > linearSolver; - resetLSIteration( ); + Eigen::Map< const Eigen::Matrix< floatType, -1, -1, Eigen::RowMajor > > J_map( getFlatJacobian( )->data( ), getNumUnknowns( ), getNumUnknowns( ) ); - while( !checkConvergence( ) && checkIteration( ) ){ + Eigen::Map< const Eigen::Vector< floatType, -1 > > R_map( getResidual( )->data( ), getNumUnknowns( ) ); - floatVector X0 = *getUnknownVector( ); + if ( *getUsePreconditioner( ) ){ - tardigradeVectorTools::solverType< floatType > linearSolver; + if( *getPreconditionerIsDiagonal( ) ){ - Eigen::Map< const Eigen::Matrix< floatType, -1, -1, Eigen::RowMajor > > J_map( getFlatJacobian( )->data( ), X0.size( ), X0.size( ) ); + Eigen::Map< const Eigen::Vector< floatType, -1 > > p_map( getFlatPreconditioner( )->data( ), getNumUnknowns( ) ); - Eigen::Map< const Eigen::Vector< floatType, -1 > > R_map( getResidual( )->data( ), X0.size( ) ); + linearSolver = tardigradeVectorTools::solverType< floatType >( p_map.asDiagonal( ) * J_map ); - if ( *getUsePreconditioner( ) ){ + dx_map = -linearSolver.solve( p_map.asDiagonal( ) * R_map ); - if( *getPreconditionerIsDiagonal( ) ){ + } + else{ - Eigen::Map< const Eigen::Vector< floatType, -1 > > p_map( getFlatPreconditioner( )->data( ), X0.size( ) ); + Eigen::Map< const Eigen::Matrix< floatType, -1, -1 > > p_map( getFlatPreconditioner( )->data( ), getNumUnknowns( ), getNumUnknowns( ) ); - linearSolver = tardigradeVectorTools::solverType< floatType >( p_map.asDiagonal( ) * J_map ); + linearSolver = tardigradeVectorTools::solverType< floatType >( p_map * J_map ); - dx_map = -linearSolver.solve( p_map.asDiagonal( ) * R_map ); + dx_map = -linearSolver.solve( p_map * R_map ); - } - else{ + } - Eigen::Map< const Eigen::Matrix< floatType, -1, -1 > > p_map( getFlatPreconditioner( )->data( ), X0.size( ), X0.size( ) ); + } + else{ - linearSolver = tardigradeVectorTools::solverType< floatType >( p_map * J_map ); + linearSolver = tardigradeVectorTools::solverType< floatType >( J_map ); + dx_map = -linearSolver.solve( R_map ); - dx_map = -linearSolver.solve( p_map * R_map ); + } - } + rank = linearSolver.rank( ); - } - else{ + if ( rank != getResidual( )->size( ) ){ - linearSolver = tardigradeVectorTools::solverType< floatType >( J_map ); - dx_map = -linearSolver.solve( R_map ); + TARDIGRADE_ERROR_TOOLS_CATCH( throw convergence_error( "The Jacobian is not full rank" ) ); - } + } - rank = linearSolver.rank( ); + } - if ( rank != getResidual( )->size( ) ){ + void hydraBase::solveNonLinearProblem( ){ + /*! + * Solve the non-linear problem + */ + + // Form the initial unknown vector + TARDIGRADE_ERROR_TOOLS_CATCH( initializeUnknownVector( ) ); - TARDIGRADE_ERROR_TOOLS_CATCH( throw convergence_error( "The Jacobian is not full rank" ) ); + floatVector deltaX( getNumUnknowns( ), 0 ); - } + resetLSIteration( ); + + while( !checkConvergence( ) && checkIteration( ) ){ + + floatVector X0 = *getUnknownVector( ); + + solveNewtonUpdate( deltaX ); updateUnknownVector( X0 + *getLambda( ) * deltaX ); diff --git a/src/cpp/tardigrade_hydra.h b/src/cpp/tardigrade_hydra.h index ad652b2..7168eeb 100644 --- a/src/cpp/tardigrade_hydra.h +++ b/src/cpp/tardigrade_hydra.h @@ -966,6 +966,8 @@ namespace tardigradeHydra{ void solveNonLinearProblem( ); + void solveNewtonUpdate( floatVector &X_tr ); + void setTolerance( const floatVector &tolerance ); void incrementIteration( ){ _iteration++; } diff --git a/src/cpp/tests/test_tardigrade_hydra.cpp b/src/cpp/tests/test_tardigrade_hydra.cpp index 19f22c3..40525a9 100644 --- a/src/cpp/tests/test_tardigrade_hydra.cpp +++ b/src/cpp/tests/test_tardigrade_hydra.cpp @@ -351,6 +351,12 @@ namespace tardigradeHydra{ } + static void solveNewtonUpdate( hydraBase &hydra, floatVector &deltaX ){ + + hydra.solveNewtonUpdate( deltaX ); + + } + static void solveNonLinearProblem( hydraBase &hydra ){ hydra.solveNonLinearProblem( ); @@ -3685,6 +3691,81 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_decomposeUnknownVector ){ } +BOOST_AUTO_TEST_CASE( test_hydraBase_solveNewtonUpdate ){ + + class hydraBaseMock : public tardigradeHydra::hydraBase { + + public: + + using tardigradeHydra::hydraBase::hydraBase; + + floatVector flatJacobian = { 1, 2, 3, 4, 5, 6, 7, 8, 2 }; + + floatVector residual = { 1, 2, 3 }; + + virtual const unsigned int getNumUnknowns( ) override{ return residual.size( ); } + + protected: + + virtual void initializeUnknownVector( ){ + + tardigradeHydra::unit_test::hydraBaseTester::set_residual( *this, residual ); + + tardigradeHydra::unit_test::hydraBaseTester::set_flatJacobian( *this, flatJacobian ); + + } + + + }; + + floatType time = 1.1; + + floatType deltaTime = 2.2; + + floatType temperature = 5.3; + + floatType previousTemperature = 23.4; + + floatVector deformationGradient = { 0.39293837, -0.42772133, -0.54629709, + 0.10262954, 0.43893794, -0.15378708, + 0.9615284 , 0.36965948, -0.0381362 }; + + floatVector previousDeformationGradient = { -0.21576496, -0.31364397, 0.45809941, + -0.12285551, -0.88064421, -0.20391149, + 0.47599081, -0.63501654, -0.64909649 }; + + floatVector previousStateVariables = { 0.53155137, 0.53182759, 0.63440096, 0.84943179, 0.72445532, + 0.61102351, 0.72244338, 0.32295891, 0.36178866, 0.22826323, + 0.29371405, 0.63097612, 0.09210494, 0.43370117, 0.43086276, + 0.4936851 , 0.42583029, 0.31226122, 0.42635131, 0.89338916, + 0.94416002, 0.50183668, 0.62395295, 0.1156184 , 0.31728548, + 0.41482621, 0.86630916, 0.25045537, 0.48303426, 0.98555979, + 0.51948512, 0.61289453, 0.12062867, 0.8263408 , 0.60306013, + 0.54506801, 0.34276383, 0.30412079 }; + + floatVector parameters = { 1, 2, 3, 4, 5 }; + + unsigned int numConfigurations = 4; + + unsigned int numNonLinearSolveStateVariables = 5; + + unsigned int dimension = 3; + + hydraBaseMock hydra( time, deltaTime, temperature, previousTemperature, deformationGradient, previousDeformationGradient, + { }, { }, + previousStateVariables, parameters, numConfigurations, numNonLinearSolveStateVariables, dimension ); + + floatVector answer = { 1./3, -2./3, 0 }; + + floatVector result( 3, 0 ); + + tardigradeHydra::unit_test::hydraBaseTester::initializeUnknownVector( hydra ); + tardigradeHydra::unit_test::hydraBaseTester::solveNewtonUpdate( hydra, result ); + + BOOST_CHECK( tardigradeVectorTools::fuzzyEquals( result, answer ) ); + +} + BOOST_AUTO_TEST_CASE( test_hydraBase_solveNonLinearProblem ){ class hydraBaseMock : public tardigradeHydra::hydraBase {