Material Point Solver for Hypo-Elastic Materials

Overview

The Material Point Solver (Hypo-Elastic) provides a standalone driver for simulating the response of a single material point using hypo-elastic constitutive models from the MAteRialMOdellingToolbox (Marmot) framework.

It is designed for:

  • Testing and verifying constitutive models,

  • Performing mixed strain/stress-controlled loading paths,

  • Recording full stress–strain histories including tangent matrices and internal variables.

The solver automatically handles time stepping, convergence checks, and output formatting.

Main Features

  • Supports mixed control (independent strain or stress control for each Voigt component)

  • Adaptive time stepping within each load step

  • Newton–Raphson iteration with configurable convergence tolerances

  • History recording for post-processing

  • CSV export for analysis or plotting

Solver Structure

The solver is organized around a few key data structures:

  • SolverOptions — defines numerical tolerances and iteration limits.

  • Step — represents a loading phase with target increments, control flags, and time control.

  • Increment — represents a sub-step automatically generated within each Step.

  • HistoryEntry — records time, stress, strain, tangent, and state variables at each increment.

Usage Workflow

A typical workflow consists of:

  1. Initialize the solver with a material model name and its properties.

  2. Set the initial state (stress and internal variables).

  3. Define one or more loading Steps with control flags and targets.

  4. Add the steps to the solver.

  5. Solve the problem.

  6. Inspect or export the resulting history.

Example Usage

#include "Marmot/MarmotMaterialPointSolverHypoElastic.h"
#include <Eigen/Dense>
#include <iostream>

int main() {
  using Vec6 = Marmot::Vector6d;

  // 1) Define the material model (must be registered in Marmot)
  std::string materialName = "LINEARELASTIC";
  double properties[] = { 210e9, 0.3 }; // Example: Young's modulus, Poisson's ratio
  int nProps = 2;

  // 2) Configure solver options
  MarmotMaterialPointSolverHypoElastic::SolverOptions options;
  options.maxIterations       = 25;
  options.residualTolerance   = 1e-10;
  options.correctionTolerance = 1e-10;

  // 3) Create solver instance
  MarmotMaterialPointSolverHypoElastic solver(materialName, properties, nProps, options);

  // 4) Set the initial state
  int nSV = 0;
  solver.getNumberOfStateVariables(nSV);
  Eigen::VectorXd initialSV = Eigen::VectorXd::Zero(nSV);
  solver.setInitialState(Vec6::Zero(), initialSV);

  // 5) Define a loading step (uniaxial strain-controlled on ε_11)
  MarmotMaterialPointSolverHypoElastic::Step step;
  step.timeStart = 0.0;
  step.timeEnd   = 1.0;
  step.dTStart   = 0.2;
  step.dTMin     = 1e-6;
  step.dTMax     = 0.5;
  step.maxIncrements = 50;

  step.strainIncrementTarget = Vec6::Zero();
  step.stressIncrementTarget = Vec6::Zero();

  step.strainIncrementTarget[0] = 1e-3; // Total applied strain over the step

  // Control flags: exactly one of each must be true
  step.isStrainComponentControlled = Eigen::Vector<bool,6>::Constant(false);
  step.isStressComponentControlled = Eigen::Vector<bool,6>::Constant(false);
  step.isStrainComponentControlled[0] = true;  // ε_11 controlled
  for (int i = 1; i < 6; ++i)
    step.isStressComponentControlled[i] = true; // others stress-controlled (zero target)

  // 6) Run the solver
  solver.addStep(step);
  solver.solve();

  // 7) Output the results
  solver.printHistory();
  solver.exportHistoryToCSV("mp_history.csv");

  return 0;
}

Typical Output

Running the above example will print information about each load step and increment, for example:

Solving step from 0 to 1
  Solving increment 1, time: 0 to 0.2, dT: 0.2
    Iteration 0, ||ddE||: 1.0e+00, ||R||: 2.5e-05
    Iteration 1, ||ddE||: 3.1e-07, ||R||: 1.2e-10
  Converged after 2 iterations.
Material Point History:
Time: 1.000000e+00
  Stress: 2.100000e+08  0.000000e+00  0.000000e+00  0.000000e+00  0.000000e+00  0.000000e+00
  Strain: 1.000000e-03  0.000000e+00  0.000000e+00  0.000000e+00  0.000000e+00  0.000000e+00
  State Vars: ...

The history can be exported as a CSV file (mp_history.csv) for visualization or analysis.

Guidelines and Notes

  • Each of the six Voigt components must be either strain-controlled or stress-controlled. Step::checkControl() is used to verify consistency when adding a step.

  • The total target increments are distributed automatically over time according to dT.

  • If convergence fails, the solver halves dT until reaching dTMin.

  • The exported CSV includes: time, stress[6], strain[6], and all state variables.

Practical Applications

  • Material verification and constitutive model calibration

  • Generating reference stress–strain curves

  • Investigating convergence and stability of material models

  • Educational demonstrations of hypo-elastic material response

class MarmotMaterialPointSolverHypoElastic

Collaboration diagram for MarmotMaterialPointSolverHypoElastic:

digraph {
    graph [bgcolor="#00000000"]
    node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
    edge [color="#1414CE"]
    "4" [label="MarmotMaterial" tooltip="MarmotMaterial"]
    "2" [label="MarmotMaterialHypoElastic" tooltip="MarmotMaterialHypoElastic"]
    "3" [label="MarmotMaterialMechanical" tooltip="MarmotMaterialMechanical"]
    "1" [label="MarmotMaterialPointSolverHypoElastic" tooltip="MarmotMaterialPointSolverHypoElastic" fillcolor="#BFBFBF"]
    "5" [label="MarmotMaterialPointSolverHypoElastic::SolverOptions" tooltip="MarmotMaterialPointSolverHypoElastic::SolverOptions"]
    "2" -> "3" [dir=forward tooltip="public-inheritance"]
    "3" -> "4" [dir=forward tooltip="public-inheritance"]
    "1" -> "2" [dir=forward tooltip="usage"]
    "1" -> "5" [dir=forward tooltip="usage"]
}

Solver for material point problems with hypo-elastic materials.

This class implements a solver for material point problems using hypo-elastic material models. It supports loading steps with controlled strain and stress components, adaptive time stepping, and history recording.

Public Functions

MarmotMaterialPointSolverHypoElastic(std::string &materialName, double *materialProperties, int nMaterialProperties, const SolverOptions &options)

Constructor for the MarmotMaterialPointSolverHypoElastic class.

Parameters:
  • materialName – Name of the hypo-elastic material model

  • materialProperties – Array of material properties

  • nMaterialProperties – Number of material properties

void addStep(const Step &step)

Add a loading step to the solver.

Parameters:

step – The Step to be added

inline std::vector<Step> getSteps() const

Get the list of added loading steps.

Returns:

A vector of Step containing the added steps

inline void clearSteps()

Clear all added loading steps.

void setInitialState(const Marmot::Vector6d &initialStress, const Eigen::VectorXd &initialStateVars)

Set the initial state of the material model.

Parameters:
  • initialStress – The initial stress in Voigt notation

  • initialStateVars – The initial state variables

inline int getNumberOfStateVariables(int &nStateVarsOut) const

Get the number of state variables in the material model.

Returns:

The number of state variables

void resetToInitialState()

Reset the solver to the initial state.

This function resets the initial stress and state variables of the material model.

void solve()

Solve the material point problem for all added steps.

inline std::vector<HistoryEntry> getHistory() const

Get the recorded history of the simulation.

Returns:

A vector of HistoryEntry containing the recorded history

inline void clearHistory()

Clear the recorded history.

void printHistory()

Print the recorded history to the console.

void exportHistoryToCSV(const std::string &filename)

Export the recorded history to a CSV file.

Parameters:

filename – The name of the CSV file to export to

Private Functions

void solveStep(const Step &step)

Solve a single loading step.

This function iterates over the increments defined in the step, calling solveIncrement for each increment. It manages time stepping and ensures that the entire step is covered.

Parameters:

step – The Step to be solved

void solveIncrement(const Increment &increment)

Solve a single increment within a loading step.

This function implements a Newton-Raphson iterative solver to compute the stress and strain state for the given increment. It updates the material state variables and records the history after convergence.

Parameters:

increment – The Increment to be solved

Throws:

std::runtime_error – if the solver does not converge

Marmot::Vector6d computeResidual(const Marmot::Vector6d &stressIncrement, const Marmot::Vector6d &target, const Increment &increment)

Compute the residual for the current increment.

This function calculates the residual vector based on the difference between the computed stress increment and the target increment, taking into account which components are controlled by strain or stress.

Parameters:
  • stressIncrement – The computed stress increment

  • target – The target (mixed) stress/strain increment

  • increment – The Increment containing control information

Returns:

The computed residual vector

void modifyTangent(Eigen::Matrix<double, 6, 6> &tangent, const Increment &increment)

Modify the material tangent matrix based on control type.

This function adjusts the tangent matrix to account for the components that are controlled by strain or stress, ensuring that the solver correctly handles mixed control scenarios. This is done by zeroing out rows corresponding to strain-controlled components and setting their diagonal entries to one.

Parameters:
  • tangent – The material tangent matrix to be modified

  • increment – The Increment containing control information

Private Members

MarmotMaterialHypoElastic *material

The hypo-elastic material model.

int nStateVars

Number of state variables in the material model.

Eigen::VectorXd stateVars

Current state variables.

Eigen::VectorXd _initialStateVars

Initial state variables.

Eigen::VectorXd stateVarsTemp

Temporary state variables for computations.

Marmot::Vector6d stress = Marmot::Vector6d::Zero()

The stress in Voigt notation.

Marmot::Vector6d _initialStress = Marmot::Vector6d::Zero()

The initial stress in Voigt notation.

Marmot::Vector6d strain = Marmot::Vector6d::Zero()

The strain in Voigt notation.

Marmot::Matrix6d dStressDStrain = Marmot::Matrix6d::Zero()

The material tangent matrix in Voigt notation.

std::vector<Step> steps

List of loading steps.

std::vector<HistoryEntry> history

History of the simulation.

const SolverOptions options

Solver options.

struct HistoryEntry

Struct to record the history of the simulation.

Each entry contains time, stress, strain, and state variables.

Public Functions

inline void print() const

Public Members

double time

Time at the history entry.

Marmot::Vector6d stress

Stress at the history entry.

Marmot::Vector6d strain

Strain at the history entry.

Marmot::Matrix6d dStressdStrain

Material tangent at the history entry.

Eigen::VectorXd stateVars

State variables at the history entry.

struct Increment

Struct to define a loading increment.

Each increment contains strain and stress increments, control flags, time information, and iteration limits.

Public Members

Marmot::Vector6d strainIncrement

Strain increment for the increment.

Marmot::Vector6d stressIncrement

Stress increment for the increment.

Eigen::Vector<bool, 6> isStrainComponentControlled

Flags to indicate which strain components are controlled.

Eigen::Vector<bool, 6> isStressComponentControlled

Flags to indicate which stress components are controlled.

double timeOld

Old time at the beginning of the increment.

double dT

Time step size for the increment.

struct SolverOptions

Struct to define solver options.

Contains parameters for controlling the solver’s behavior.

Public Members

int maxIterations = 25

Maximum number of iterations per increment.

double residualTolerance = 1e-10

Convergence tolerance.

double correctionTolerance = 1e-10

Correction tolerance.

struct Step

Struct to define a loading step.

Each step contains target strain and stress states, time information and time step control parameters.

Public Functions

inline void checkControl() const

Check that for each component, either strain or stress is controlled.

Throws:

std::runtime_error – if the condition is not met

Public Members

Marmot::Vector6d strainIncrementTarget

Target strain increment for the step.

Marmot::Vector6d stressIncrementTarget

Target stress increment for the step.

Eigen::Vector<bool, 6> isStrainComponentControlled

Flags to indicate which strain components are controlled.

Eigen::Vector<bool, 6> isStressComponentControlled

Flags to indicate which stress components are controlled.

double timeStart = 0.0

Start time of the step.

double timeEnd = 1.0

End time of the step.

double dTStart = 0.1

Initial time step size.

double dTMin = 1e-6

Minimum time step size.

double dTMax = 0.5

Maximum time step size.

int maxIncrements = 100

Maximum number of increments in the step.