|
NEML2 2.0.0
|
One of the most notable differences between constitutive models and feed-forward neural networks is that updating certain stiff systems with implicit methods is often more computationally efficient compared to explicit algorithms.
A generally nonlinear, recursive, implicit system of equations take the following form
\begin{align*} \mathbf{r}(\tilde{\mathbf{s}}) & = f(\tilde{\mathbf{s}}, \mathbf{f}, \mathbf{s}_n, \mathbf{f}_n; \mathbf{p}), \\ \mathbf{s} &= \mathop{\mathrm{root}}\limits_{\tilde{\mathbf{s}}} (\mathbf{r}). \end{align*}
Here \( \mathbf{r} \) represents the residual (for root-fiding) of the system of equations, and \( \mathbf{s} \), \( \mathbf{f} \), \( \mathbf{s}_n \), \( \mathbf{f}_n \) are defined by four reserved sub-axes in NEML2:
state sub-axis hosts input variables in set \( \tilde{\mathbf{s}} \). The variables on this sub-axis are the primary unknowns to be solved for. After solving the system, the state sub-axis hosts output variables in set \( \mathbf{s} \).forces sub-axis hosts prescribed input variables in set \( \mathbf{f} \). These variables are prescribed and, by definition, do not change while the system is being solved.old_state and old_forces sub-axes respectively correspond to \( \mathbf{s}_n \) and \( \mathbf{f}_n \). These variables correspond to the previous solution to the system to facilitate the recursive definition of internal variables in history-dependent models. The equivalent plastic strain in plasticity models is a well-known example.The Perzyna viscoplasticity model mentioned in a previous tutorial takes this form:
\begin{align} \boldsymbol{\varepsilon}^e & = \boldsymbol{\varepsilon} - \boldsymbol{\varepsilon}^p, \label{1} \\ \boldsymbol{\sigma} & = 3K\operatorname{vol}\boldsymbol{\varepsilon}^e + 2G\operatorname{dev}\boldsymbol{\varepsilon}^e, \label{2} \\ \bar{\sigma} & = J_2(\boldsymbol{\sigma}), \label{3} \\ f^p & = \bar{\sigma} - \sigma_y, \label{4} \\ \boldsymbol{N} & = \pdv{f^p}{\boldsymbol{\sigma}}, \label{5} \\ \dot{\gamma} & = \left( \dfrac{\left< f^p \right>}{\eta} \right)^n, \label{6} \\ \dot{\boldsymbol{\varepsilon}}^p & = \dot{\gamma} \boldsymbol{N}. \label{7} \end{align}
Its residual can be defined using backward-Euler time integration as
\begin{align} \mathbf{r}\left( \boldsymbol{\varepsilon}^p \right) = \boldsymbol{\varepsilon}^p - \boldsymbol{\varepsilon}^p_n - \left( t - t_n \right) \dot{\boldsymbol{\varepsilon}}^p, \label{8} \end{align}
the solution procedure of which is often referred to as the return-mapping algorithm in the solid mechanics community.
This tutorial illustrates the use of ImplicitUpdate in conjunction with a Newton-Raphson solver to perform the implicit update.
After searching among existing models provided by NEML2, we are able to translate each of these equations into a NEML2 model. The input file looks like
Note the one-to-one correspondance between the models and the equations.
The structure of the system of equations can be summarized using the code below.
C++
Output: @list-output:ex1
Python
Output: @list-output:ex2
Once the system of equations are properly defined, we can use the ImplicitUpdate to solve the system of equations. The ImplicitUpdate is responsible for the following:
NEML2 offers three fully vectorized Newton solvers to be used in conjunction with ImplicitUpdate:
In addition, the assembly routines as well as the application of the implicit function theorem are also implemented in a vectorized fashion.
The additional sections needed in the input file are
The ImplicitUpdate model can then be invoked in the same way as regular models.
C++
Output: @list-output:ex3
Python
Output: @list-output:ex4
Unlike other regular models, declaring variables on the correct sub-axes is important because NEML2 relies on the reserved sub-axes (state, forces, etc.)
forces sub-axis are skipped because they are not required in the assembly of the linear system.As long as models are defined using the correct sub-axis definitions and satisfy some mild continuity requirements, NEML2 guarantees the correctness of the variable derivatives after one or more implicit updates, up to machine precision. The same guarantee also applies to user-defined custom models.
This is a significant advantage compared to some of the alternative constitutive model libraries, especially in the context of coupling with external PDE solvers. For example, in the context of finite element method, thermodynamic forces (e.g. strain) are calculated at each quadrature point, and the constitutive library (e.g. NEML2) is responsible for updating the thermodynamic state variables (e.g. stress, plastic strain, etc.), which are then used in the residual definition of the discretized PDE. Therefore, the exact derivatives of the state variables with respect to the forces are the key to the assembly of the exact Jacobian of the descretized PDE, which is in turn the fundamental requirement for optimal convergence for many nonlinear solvers.