diff --git a/docs/sphinx/changelog.rst b/docs/sphinx/changelog.rst index 5e19dbc..7bc07d4 100644 --- a/docs/sphinx/changelog.rst +++ b/docs/sphinx/changelog.rst @@ -20,6 +20,7 @@ New Features - Added a relaxed solver as the default (:pull:`172`). By `Nathan Miller`_. - Added ability to turn on additional messages for failed solves (:pull:`177`). By `Nathan Miller`_. - Added the ability to scale the incoming load information by a scale factor (:pull:`178`). By `Nathan Miller`_. +- Added the ability to sub-cycle the analysis to try and improve convergence (:pull:`179`). By `Nathan Miller`_. Internal Changes ================ diff --git a/src/cpp/tardigrade_hydra.cpp b/src/cpp/tardigrade_hydra.cpp index ff69589..a3fa05e 100644 --- a/src/cpp/tardigrade_hydra.cpp +++ b/src/cpp/tardigrade_hydra.cpp @@ -44,6 +44,7 @@ namespace tardigradeHydra{ const unsigned int maxLSIterations, const floatType lsAlpha, const bool use_preconditioner, const unsigned int preconditioner_type ) : _dimension( dimension ), _configuration_unknown_count( configuration_unknown_count ), + _stress_size( configuration_unknown_count ), _time( time ), _deltaTime( deltaTime ), _temperature( temperature ), _previousTemperature( previousTemperature ), _deformationGradient( deformationGradient ), @@ -224,22 +225,19 @@ namespace tardigradeHydra{ } - void hydraBase::decomposeUnknownVector( ){ + void hydraBase::updateConfigurationsFromUnknownVector( ){ /*! - * Decompose the unknown vector into the cauchy stress, configurations, and state variables used for the non-linear solve + * Update the configurations from the unknown vector */ const floatVector *unknownVector = getUnknownVector( ); - // Set the stress - extractStress( ); - // Set the configurations auto configurations = get_setDataStorage_configurations( ); auto inverseConfigurations = get_setDataStorage_inverseConfigurations( ); - computeConfigurations( unknownVector, getStress( )->size( ), *getDeformationGradient( ), *configurations.value, *inverseConfigurations.value ); + computeConfigurations( unknownVector, *getStressSize( ), *getDeformationGradient( ), *configurations.value, *inverseConfigurations.value ); // Extract the remaining state variables required for the non-linear solve auto nonLinearSolveStateVariables = get_setDataStorage_nonLinearSolveStateVariables( ); @@ -254,6 +252,18 @@ namespace tardigradeHydra{ } + void hydraBase::decomposeUnknownVector( ){ + /*! + * Decompose the unknown vector into the cauchy stress, configurations, and state variables used for the non-linear solve + */ + + // Set the stress + extractStress( ); + + updateConfigurationsFromUnknownVector( ); + + } + void hydraBase::decomposeStateVariableVector( ){ /*! * Decompose the incoming state variable vector setting the different configurations along the way @@ -2034,9 +2044,177 @@ namespace tardigradeHydra{ } - void hydraBase::evaluate( ){ + void hydraBase::setScaledQuantities( ){ + /*! + * Set the scaled quantities + */ + + _scaled_time = ( _scale_factor - 1 ) * _deltaTime + _time; + + _scaled_deltaTime = _scale_factor * _deltaTime; + + _scaled_temperature = _scale_factor * ( _temperature - _previousTemperature ) + _previousTemperature; + + _scaled_deformationGradient = _scale_factor * ( _deformationGradient - _previousDeformationGradient ) + _previousDeformationGradient; + + _scaled_additionalDOF = _scale_factor * ( _additionalDOF - _previousAdditionalDOF ) + _previousAdditionalDOF; + + } + + void hydraBase::setScaleFactor( const floatType &value ){ + /*! + * Set the value of the scale factor. Will automatically re-calculate the deformation and trial stresses + * + * \param &value: The value of the scale factor + */ + + // Update the scale factor + _scale_factor = value; + + // Update the scaled quantities + setScaledQuantities( ); + + // Copy the current unknown vector + floatVector unknownVector = *getUnknownVector( ); + + // Reset the iteration data + resetIterationData( ); + + // Set the unknown vector + setX( unknownVector ); + + // Update the deformation quantities + updateConfigurationsFromUnknownVector( ); + + // Compute the new trial stress + std::copy( getStress( )->begin( ), getStress( )->end( ), unknownVector.begin( ) ); + + // Re-set the unknown vector + setX( unknownVector ); + + // Extract the stress + extractStress( ); + + } + + const bool hydraBase::allowStepGrowth( const unsigned int &num_good ){ + /*! + * Function to determine if we can increase the step-size for the sub-cycler + * + * \param &num_good: The number of good increments since the last failure + */ + + if ( num_good >= ( *getNumGoodControl( ) ) ){ + + return true; + + } + + return false; + + } + + void hydraBase::evaluate( const bool &use_subcycler ){ + /*! + * Solver the non-linear problem and update the variables + * + * \param &use_subcycler: Flag for if the subcycler should be used for difficult analyses (defaults to false) + */ + + try{ + + evaluateInternal( ); + + return; + + } + catch( std::exception &e ){ + + if ( !use_subcycler ){ + + throw; + + } + + if ( ( *getFailureVerbosityLevel( ) ) > 0 ){ + addToFailureOutput( "\n\n" ); + addToFailureOutput( "#########################################\n" ); + addToFailureOutput( "### ENTERING SUB-CYCLER ###\n" ); + addToFailureOutput( "#########################################\n" ); + addToFailureOutput( "\n\n" ); + } + + floatType sp = 0.0; + + floatType ds = ( *getCutbackFactor( ) ); + + unsigned int num_good = 0; + + // Set the unknown vector to the initial unknown. We're using setX because we call setScaleFactor right away which will update the unknown vector + setX( _initialX ); + + while ( sp < 1.0 ){ + + try{ + + if ( ( *getFailureVerbosityLevel( ) ) > 0 ){ + addToFailureOutput( "\n\n" ); + addToFailureOutput( "######### PSEUDO-TIME INCREMENT #########\n" ); + addToFailureOutput( "\n\n sp, ds: " + std::to_string( sp ) + ", " + std::to_string( ds ) ); + addToFailureOutput( "\n" ); + } + + setScaleFactor( sp + ds ); // Update the scaling factor + + resetIterations( ); // Reset the non-linear iteration count + + evaluateInternal( ); // Try to solve the non-linear problem + + sp += ds; // Update the pseudo-time + + num_good++; // Update the number of good iterations + + // Grow the step if possible + if ( allowStepGrowth( num_good ) ){ + + ds *= ( *getGrowthFactor( ) ); + + } + + // Make sure s will be less than or equal to 1 + if ( sp + ds > 1.0 ){ + + ds = 1.0 - sp; + + } + + } + catch( std::exception &e ){ + + // Reduce the time-step and try again + num_good = 0; + + ds *= ( *getCutbackFactor( ) ); + + setX( _initialX ); // Reset X to the last good point + + if ( ds < *getMinDS( ) ){ + + throw; + + } + + } + + } + + } + + } + + void hydraBase::evaluateInternal( ){ /*! - * Solve the non-linear problem and update the variables + * Solve the non-linear problem with the current scaling and update the variables */ // Reset the counters for the number of steps being performed diff --git a/src/cpp/tardigrade_hydra.h b/src/cpp/tardigrade_hydra.h index db289e3..f653500 100644 --- a/src/cpp/tardigrade_hydra.h +++ b/src/cpp/tardigrade_hydra.h @@ -1080,6 +1080,9 @@ namespace tardigradeHydra{ //! Get a reference to the number of unknowns in each configuration const unsigned int* getConfigurationUnknownCount( ){ return &_configuration_unknown_count; } + //! Get a reference to the number of components of the stress + const unsigned int* getStressSize( ){ return &_stress_size; } + //! Get a reference to the current time const floatType* getTime( ){ return getScaledTime( ); } @@ -1413,7 +1416,7 @@ namespace tardigradeHydra{ const floatVector* getPreviousStress( ); - virtual void evaluate( ); + virtual void evaluate( const bool &use_subcycler = false ); virtual void computeTangents( ); @@ -1494,8 +1497,8 @@ namespace tardigradeHydra{ //! Get a scale factor for the deformation const floatType *getScaleFactor( ){ return &_scale_factor; } - //! Set the value of the scale factor. Will automatically re-calculate the deformation - const void setScaleFactor( const floatType &value ){ _scale_factor = value; setScaledQuantities( ); updateUnknownVector( *getUnknownVector( ) ); } + //! Set the value of the scale factor. Will automatically re-calculate the deformation and trial stresses + void setScaleFactor( const floatType &value ); const floatType *getScaledTime( ){ return &_scaled_time; } @@ -1507,6 +1510,24 @@ namespace tardigradeHydra{ const floatVector *getScaledAdditionalDOF( ){ return &_scaled_additionalDOF; } + const floatType *getCutbackFactor( ){ return &_cutback_factor; } + + const unsigned int *getNumGoodControl( ){ return &_num_good_control; } + + const floatType *getGrowthFactor( ){ return &_growth_factor; } + + const floatType *getMinDS( ){ return &_minDS; } + + void setCutbackFactor( const floatType &value ){ _cutback_factor = value; } + + void setNumGoodControl( const unsigned int &value ){ _num_good_control = value; } + + void setGrowthFactor( const floatType &value ){ _growth_factor = value; } + + void setMinDS( const floatType &value ){ _minDS = value; } + + const bool allowStepGrowth( const unsigned int &num_good ); + protected: // Setters that the user may need to access but not override @@ -1521,6 +1542,8 @@ namespace tardigradeHydra{ virtual void extractStress( ); + virtual void updateConfigurationsFromUnknownVector( ); + virtual void decomposeUnknownVector( ); virtual void decomposeStateVariableVector( ); @@ -1541,6 +1564,8 @@ namespace tardigradeHydra{ virtual void initializePreconditioner( ); + virtual void evaluateInternal( ); + //! Update the line-search lambda parameter virtual void updateLambda( ){ _lambda *= 0.5; } @@ -1553,19 +1578,7 @@ namespace tardigradeHydra{ virtual void performGradientStep( const floatVector &X0 ); //! Update the scaled quantities - virtual void setScaledQuantities( ){ - - _scaled_time = ( _scale_factor - 1 ) * _deltaTime + _time; - - _scaled_deltaTime = _scale_factor * _deltaTime; - - _scaled_temperature = _scale_factor * ( _temperature - _previousTemperature ) + _previousTemperature; - - _scaled_deformationGradient = _scale_factor * ( _deformationGradient - _previousDeformationGradient ) + _previousDeformationGradient; - - _scaled_additionalDOF = _scale_factor * ( _additionalDOF - _previousAdditionalDOF ) + _previousAdditionalDOF; - - } + virtual void setScaledQuantities( ); const floatType *get_baseResidualNorm( ); @@ -1776,6 +1789,8 @@ namespace tardigradeHydra{ unsigned int _configuration_unknown_count; //!< The number of unknowns required for a configuration. Used to ensure that the unknown and state variable vectors are the right size. Must be set by all inheriting classes. For 3D classical continuum this will be 9, for higher order theories this will change. + unsigned int _stress_size; //!< The number of terms in the stress measures. For 3D classical continuum this will be 9, for higher order theories this will change. + floatType _time; //!< The current time floatType _deltaTime; //!< The change in time @@ -1948,6 +1963,14 @@ namespace tardigradeHydra{ floatType _scale_factor = 1.0; //!< A scale factor applied to the incoming loading (deformation, temperature, etc.) + floatType _cutback_factor = 0.5; //!< The factor by which the pseudo-time will be scaled if a solve fails + + floatType _growth_factor = 1.2; //!< The factor by which the pseudo-time will be scaled if we can grow the pseudo-timestep + + unsigned int _num_good_control = 2; //!< The number of good iterations we need to have before we try and increase the timestep + + floatType _minDS = 1e-2; //!< The minimum allowable pseudo-timestep + TARDIGRADE_HYDRA_DECLARE_ITERATION_STORAGE( private, configurations, floatVector, passThrough ) TARDIGRADE_HYDRA_DECLARE_PREVIOUS_STORAGE( private, previousConfigurations, floatVector, passThrough ) diff --git a/src/cpp/tardigrade_hydraMicromorphic.cpp b/src/cpp/tardigrade_hydraMicromorphic.cpp index ab60800..4f2cc50 100644 --- a/src/cpp/tardigrade_hydraMicromorphic.cpp +++ b/src/cpp/tardigrade_hydraMicromorphic.cpp @@ -63,6 +63,19 @@ namespace tardigradeHydra{ } + void hydraBaseMicromorphic::setScaledQuantities( ){ + /*! + * Scale the current values by the scale factor + */ + + hydraBase::hydraBase::setScaledQuantities( ); + + _scaled_microDeformation = ( *getScaleFactor( ) ) * ( _microDeformation - _previousMicroDeformation ) + _previousMicroDeformation; + + _scaled_gradientMicroDeformation = ( *getScaleFactor( ) ) * ( _gradientMicroDeformation - _previousGradientMicroDeformation ) + _previousGradientMicroDeformation; + + } + void hydraBaseMicromorphic::initializeUnknownVector( ){ /*! * Initialize the unknown vector for the non-linear solve. @@ -179,39 +192,50 @@ namespace tardigradeHydra{ } - void hydraBaseMicromorphic::decomposeUnknownVector( ){ + void hydraBaseMicromorphic::updateConfigurationsFromUnknownVector( ){ /*! * Decompose the incoming unknown vector setting the different configurations along the way - * - * The state variable vector is assumed to be of the form: - * - * \f$ \text{ISV} = \left\{\bf{S}^2, \bf{\Sigma}, \bf{M}, \bf{F}^2, \bf{F}^3, \cdots, \bf{F}^n, \bf{\chi}^2, \bf{\chi}^3, \cdots, \bf{\chi}^n \frac{\partial}{\partial \bf{X}} \bf{\chi}^2, \frac{\partial}{\partial \bf{X}} \bf{\chi}^3, \cdots, \frac{\partial}{\partial \bf{X}} \bf{\chi}^n, \xi^1, \xi^2, \cdots, \xi^m\right\} \f$ - * - * where \f$\bf{S}^2\f$ is the second Piola Kirchhoff stress, \f$\bf{\Sigma}\f$ is the reference symmetric micro stress, and \f$\bf{M}\f$ - * is the reference higher order stress the \f$\bf{F}\f$ are the different deformation gradients (configurations), \f$\bf{\chi}\f$ are the micro-deformations, - * \f$\xi^y\f$ are the other variables to be solved during the non-linear solve, and \f$\eta^z\f$ are other state variables. Note - * that we decompose the deformation gradient and micro-deformation as - * - * \f$\bf{F} = \bf{F}^1 \bf{F}^2 \cdots \bf{F}^n\f$ - * - * \f$\bf{\chi} = \bf{\chi}^1 \bf{\chi}^2 \cdots \bf{\chi}^n\f$ - * - * and so because \f$\bf{F}\f$ and \f$\bf{\chi}\f$ are provided we can solve for \f$\bf{F}^1\f$ and \f$\bf{\chi}\f$. Typically, - * this configuration would be the elastic configuration (i.e., the configuration that generates the stress) though we do not insist that users follow convention. - * - * NOTE: Though we overload the decomposeStateVariableVector in hydraBase this function will not be called in hydraBase's constructor because - * virtual functions do not come into being during the construction of parent classes constructors. We could work around this but instead - * we will overload and define a local method to do the decomposition of the micro-deformation tensors. */ - // Call the parent class decomposition - hydraBase::decomposeUnknownVector( ); + hydraBase::updateConfigurationsFromUnknownVector( ); - // Decompose the micro-deformation decomposeUnknownVectorMicroConfigurations( ); } +// void hydraBaseMicromorphic::decomposeUnknownVector( ){ +// /*! +// * Decompose the incoming unknown vector setting the different configurations along the way +// * +// * The state variable vector is assumed to be of the form: +// * +// * \f$ \text{ISV} = \left\{\bf{S}^2, \bf{\Sigma}, \bf{M}, \bf{F}^2, \bf{F}^3, \cdots, \bf{F}^n, \bf{\chi}^2, \bf{\chi}^3, \cdots, \bf{\chi}^n \frac{\partial}{\partial \bf{X}} \bf{\chi}^2, \frac{\partial}{\partial \bf{X}} \bf{\chi}^3, \cdots, \frac{\partial}{\partial \bf{X}} \bf{\chi}^n, \xi^1, \xi^2, \cdots, \xi^m\right\} \f$ +// * +// * where \f$\bf{S}^2\f$ is the second Piola Kirchhoff stress, \f$\bf{\Sigma}\f$ is the reference symmetric micro stress, and \f$\bf{M}\f$ +// * is the reference higher order stress the \f$\bf{F}\f$ are the different deformation gradients (configurations), \f$\bf{\chi}\f$ are the micro-deformations, +// * \f$\xi^y\f$ are the other variables to be solved during the non-linear solve, and \f$\eta^z\f$ are other state variables. Note +// * that we decompose the deformation gradient and micro-deformation as +// * +// * \f$\bf{F} = \bf{F}^1 \bf{F}^2 \cdots \bf{F}^n\f$ +// * +// * \f$\bf{\chi} = \bf{\chi}^1 \bf{\chi}^2 \cdots \bf{\chi}^n\f$ +// * +// * and so because \f$\bf{F}\f$ and \f$\bf{\chi}\f$ are provided we can solve for \f$\bf{F}^1\f$ and \f$\bf{\chi}\f$. Typically, +// * this configuration would be the elastic configuration (i.e., the configuration that generates the stress) though we do not insist that users follow convention. +// * +// * NOTE: Though we overload the decomposeStateVariableVector in hydraBase this function will not be called in hydraBase's constructor because +// * virtual functions do not come into being during the construction of parent classes constructors. We could work around this but instead +// * we will overload and define a local method to do the decomposition of the micro-deformation tensors. +// */ +// +// // Call the parent class decomposition +// hydraBase::decomposeUnknownVector( ); +// +// // Decompose the micro-deformation +// decomposeUnknownVectorMicroConfigurations( ); +// +// } + void hydraBaseMicromorphic::computeGradientMicroConfigurations( const floatVector *data_vector, unsigned int start_index, const floatVector &configurations, const floatVector µConfigurations, const floatVector &gradientMicroConfiguration, floatVector &gradientMicroConfigurations ){ @@ -245,7 +269,7 @@ namespace tardigradeHydra{ const unsigned int sot_dim = getSOTDimension( ); const unsigned int num_configs = *getNumConfigurations( ); - unsigned int start_index = getStress( )->size( ) + ( num_configs - 1 ) * sot_dim; + unsigned int start_index = ( *getStressSize( ) ) + ( num_configs - 1 ) * sot_dim; auto microConfigurations = get_setDataStorage_microConfigurations( ); diff --git a/src/cpp/tardigrade_hydraMicromorphic.h b/src/cpp/tardigrade_hydraMicromorphic.h index aec0c53..687a2da 100644 --- a/src/cpp/tardigrade_hydraMicromorphic.h +++ b/src/cpp/tardigrade_hydraMicromorphic.h @@ -97,7 +97,9 @@ namespace tardigradeHydra{ //Utility functions virtual void initializeUnknownVector( ) override; - virtual void decomposeUnknownVector( ) override; + virtual void updateConfigurationsFromUnknownVector( ) override; + +// virtual void decomposeUnknownVector( ) override; virtual void decomposeUnknownVectorMicroConfigurations( ); @@ -105,18 +107,7 @@ namespace tardigradeHydra{ virtual void decomposeStateVariableVectorMicroConfigurations( ); - virtual void setScaledQuantities( ) override{ - /*! - * Scale the current values by the scale factor - */ - - hydraBase::hydraBase::setScaledQuantities( ); - - _scaled_microDeformation = ( *getScaleFactor( ) ) * ( _microDeformation - _previousMicroDeformation ) + _previousMicroDeformation; - - _scaled_gradientMicroDeformation = ( *getScaleFactor( ) ) * ( _gradientMicroDeformation - _previousGradientMicroDeformation ) + _previousGradientMicroDeformation; - - } + virtual void setScaledQuantities( ) override; private: diff --git a/src/cpp/tests/test_tardigrade_hydra.cpp b/src/cpp/tests/test_tardigrade_hydra.cpp index bd1a00f..d33efdf 100644 --- a/src/cpp/tests/test_tardigrade_hydra.cpp +++ b/src/cpp/tests/test_tardigrade_hydra.cpp @@ -6085,7 +6085,7 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_updateUnknownVector, * boost::unit_test::to } -BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ +BOOST_AUTO_TEST_CASE( test_hydraBase_evaluateInternal, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ class residualMock : public tardigradeHydra::residualBase{ @@ -6118,6 +6118,8 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate, * boost::unit_test::tolerance( DE void setInitialX( ){ _initialX = _mockInitialX; } + void public_evaluateInternal( ){ evaluateInternal( ); } + protected: floatVector _mockInitialX = { 1, 1, 1, 1, 1, 1, 1, 1, 1, @@ -6171,7 +6173,7 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate, * boost::unit_test::tolerance( DE hydra.setUseRelaxedSolve( false ); - BOOST_CHECK_THROW( hydra.evaluate( ), tardigradeHydra::convergence_error ); + BOOST_CHECK_THROW( hydra.public_evaluateInternal( ), tardigradeHydra::convergence_error ); BOOST_TEST( ( *hydra.getUseLevenbergMarquardt( ) ) ); @@ -6181,7 +6183,7 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate, * boost::unit_test::tolerance( DE } -BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate2, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ +BOOST_AUTO_TEST_CASE( test_hydraBase_evaluateInternal2, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ class residualMock : public tardigradeHydra::residualBase{ @@ -6216,6 +6218,8 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate2, * boost::unit_test::tolerance( D void setInitialX( ){ _initialX = _mockInitialX; } + void public_evaluateInternal( ){ evaluateInternal( ); } + protected: floatVector _mockInitialX = { 1, 1, 1, 1, 1, 1, 1, 1, 1, @@ -6273,7 +6277,7 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate2, * boost::unit_test::tolerance( D { }, { }, previousStateVariables, parameters, numConfigurations, numNonLinearSolveStateVariables, dimension ); - hydra.evaluate( ); + hydra.public_evaluateInternal( ); BOOST_TEST( !( *hydra.getUseLevenbergMarquardt( ) ) ); @@ -6285,6 +6289,258 @@ BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate2, * boost::unit_test::tolerance( D } +BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ + + class hydraBaseMock : public tardigradeHydra::hydraBase{ + + public: + + using tardigradeHydra::hydraBase::hydraBase; + + unsigned int num_evaluateInternalCalls = 0; + + unsigned int num_updateUnknownVectorCalls = 0; + + floatVector X = { 1, 2, 3 }; + + protected: + + virtual void updateUnknownVector( const floatVector &newX ) override{ + num_updateUnknownVectorCalls++; + } + + virtual void evaluateInternal( ) override{ + num_evaluateInternalCalls++; + } + + }; + + floatType time = 1.1; + + floatType deltaTime = 2.2; + + floatType temperature = 5.3; + + floatType previousTemperature = 23.4; + + floatVector deformationGradient = { 1.05, 0, 0, + 0.00, 1, 0, + 0.00, 1, 1}; + + floatVector previousDeformationGradient = { -0.21576496, -0.31364397, 0.45809941, + -0.12285551, -0.88064421, -0.20391149, + 0.47599081, -0.63501654, -0.64909649 }; + + floatVector previousStateVariables = { 0, 0, 0, 0, 0, 0, 0, 0, 0 }; + + floatVector parameters = { 123.4, 56.7 }; + + unsigned int numConfigurations = 2; + + unsigned int numNonLinearSolveStateVariables = 0; + + unsigned int dimension = 3; + + hydraBaseMock hydra( time, deltaTime, temperature, previousTemperature, deformationGradient, previousDeformationGradient, + { }, { }, + previousStateVariables, parameters, numConfigurations, numNonLinearSolveStateVariables, dimension ); + + hydra.evaluate( ); + + BOOST_TEST( hydra.num_evaluateInternalCalls == 1 ); + +} + +BOOST_AUTO_TEST_CASE( test_hydraBase_evaluate2, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ + + class residualBaseMock : public tardigradeHydra::residualBase{ + + using tardigradeHydra::residualBase::residualBase; + + }; + + class residualBaseMockStress : public tardigradeHydra::residualBase{ + + using tardigradeHydra::residualBase::residualBase; + + using tardigradeHydra::residualBase::setStress; + + floatVector cauchyStress = { 1, 2, 3, 4, 5, 6, 7, 8, 9 }; + + virtual void setStress( ) override{ + + setStress( cauchyStress + ( *hydra->getDeformationGradient( ) ) ); + + } + + }; + + class hydraBaseMock : public tardigradeHydra::hydraBase{ + + public: + + using tardigradeHydra::hydraBase::hydraBase; + + using tardigradeHydra::hydraBase::setResidualClasses; + + unsigned int num_evaluateInternalCalls = 0; + + unsigned int num_updateUnknownVectorCalls = 0; + + unsigned int num_updateConfigurationsFromUnknownVectorCalls = 0; + + std::vector< unsigned int > fail_indices = { 0, 1, 3 }; + + residualBaseMockStress r1; + + residualBaseMock r2; + + residualBaseMock r3; + + unsigned int s1 = 9; + + unsigned int s2 = 10; + + unsigned int s3 = 3; + + floatVector expected_scale_factors = { 1.0, 0.5, 0.25, 0.5, 0.375, 0.50, 0.65, 0.83, 1.0 }; + + floatMatrix expected_unknownVectors = { { 2.05, 2. , 3. , 4. , 6. , 6. , 7. , 9. , 10. , + 1. , 0. , 0. , 0. , 1. , 0. , 0. , 0. , 1. , + 0.1 , 0.2 , 0.3 , 0.4 }, + { 1.41711752, 1.84317801, 3.2290497 , 3.93857225, 5.0596779 , + 5.89804426, 7.23799541, 8.18249173, 9.17545175, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.10067628, 1.76476702, 3.34357456, 3.90785837, 4.58951684, + 5.84706638, 7.35699311, 7.7737376 , 8.76317763, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.41711752, 1.84317801, 3.2290497 , 3.93857225, 5.0596779 , + 5.89804426, 7.23799541, 8.18249173, 9.17545175, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.2588969 , 1.80397252, 3.28631213, 3.92321531, 4.82459737, + 5.87255532, 7.29749426, 7.97811466, 8.96931469, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.41711752, 1.84317801, 3.2290497 , 3.93857225, 5.0596779 , + 5.89804426, 7.23799541, 8.18249173, 9.17545175, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.60698226, 1.89022461, 3.16033479, 3.95700057, 5.34177453, + 5.92863098, 7.16659678, 8.42774421, 9.42281623, 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 1.83481996, 1.94668053, 3.0778769 , 3.97911456, 5.68029048, + 5.96533505, 7.08091844, 8.72204719, 9.7196536 , 1. , + 0. , 0. , 0. , 1. , 0. , + 0. , 0. , 1. , 0.1 , 0.2 , + 0.3 , 0.4 }, + { 2.05, 2. , 3. , 4. , 6. , 6. , 7. , 9. , 10. , + 1. , 0. , 0. , 0. , 1. , 0. , 0. , 0. , 1. , + 0.1 , 0.2 , 0.3 , 0.4 } }; + + protected: + + virtual void updateUnknownVector( const floatVector &newX ) override{ + tardigradeHydra::hydraBase::updateUnknownVector( newX ); + num_updateUnknownVectorCalls++; + } + + virtual void updateConfigurationsFromUnknownVector( ) override{ + tardigradeHydra::hydraBase::updateConfigurationsFromUnknownVector( ); + num_updateConfigurationsFromUnknownVectorCalls++; + } + + virtual void evaluateInternal( ) override{ + + _initialX = *getUnknownVector( ); + + BOOST_TEST( ( *getScaleFactor( ) ) == expected_scale_factors[ num_evaluateInternalCalls ] ); + + BOOST_TEST( ( *getUnknownVector( ) ) == expected_unknownVectors[ num_evaluateInternalCalls ], CHECK_PER_ELEMENT ); + + if ( std::find( fail_indices.begin( ), fail_indices.end( ), num_evaluateInternalCalls ) != fail_indices.end( ) ){ + + num_evaluateInternalCalls++; + + throw std::runtime_error("failure in evaluateInternal"); + + } + + num_evaluateInternalCalls++; + + } + + virtual void setResidualClasses( ){ + + r1 = residualBaseMockStress( this, s1 ); + + r2 = residualBaseMock( this, s2 ); + + r3 = residualBaseMock( this, s3 ); + + std::vector< tardigradeHydra::residualBase* > residuals( 3 ); + + residuals[ 0 ] = &r1; + + residuals[ 1 ] = &r2; + + residuals[ 2 ] = &r3; + + setResidualClasses( residuals ); + + } + + }; + + floatType time = 1.1; + + floatType deltaTime = 2.2; + + floatType temperature = 5.3; + + floatType previousTemperature = 23.4; + + floatVector deformationGradient = { 1.05, 0, 0, + 0.00, 1, 0, + 0.00, 1, 1}; + + floatVector previousDeformationGradient = { -0.21576496, -0.31364397, 0.45809941, + -0.12285551, -0.88064421, -0.20391149, + 0.47599081, -0.63501654, -0.64909649 }; + + floatVector previousStateVariables = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0.1, 0.2, 0.3, 0.4 }; + + floatVector parameters = { 123.4, 56.7 }; + + unsigned int numConfigurations = 2; + + unsigned int numNonLinearSolveStateVariables = 4; + + unsigned int dimension = 3; + + hydraBaseMock hydra( time, deltaTime, temperature, previousTemperature, deformationGradient, previousDeformationGradient, + { }, { }, + previousStateVariables, parameters, numConfigurations, numNonLinearSolveStateVariables, dimension ); + + hydra.evaluate( true ); + + BOOST_TEST( hydra.num_evaluateInternalCalls == 9 ); + + BOOST_TEST( hydra.num_updateUnknownVectorCalls == 0 ); + + BOOST_TEST( hydra.num_updateConfigurationsFromUnknownVectorCalls == 8 ); + +} + BOOST_AUTO_TEST_CASE( test_setDataStorageBase, * boost::unit_test::tolerance( DEFAULT_TEST_TOLERANCE ) ){ class residualMock : public tardigradeHydra::residualBase{