LLNL / ExaCMech

BSD 3-Clause "New" or "Revised" License
17 stars 0 forks source link

Add a wrapper for Abaqus #23

Open vigneshwaran-radhakrishnan opened 4 days ago

vigneshwaran-radhakrishnan commented 4 days ago

I'm looking for a wrapper for calling ExaCMech from Abaqus/UMAT for two reasons. First, the wrapper will make ExaCMech more accessible to many users of UMAT. Nowadays, UMAT can be called from several FE/FFT codes, including ExaConsist. Second, ExaConsist is still being improved to support periodic boundary conditions, controlled triaxiality, and Lode loading conditions. While these features are being implemented, a wrapper for ExaCMech will allow users to run jobs in Abaqus using ExaCMech's material models.

rcarson3 commented 4 days ago

@vigneshwaran-radhakrishnan given I don't really use Abaqus and I last looked at the UMAT interface ~5 years ago when I added UMAT support into ExaConstit: https://github.com/LLNL/ExaConstit/blob/exaconstit-dev/src/mechanics_umat.cpp#L306-L578 . You will need to help me in creating a wrapper over ExaCMech given my lack of expertise around the Abaqus / UMAT interfaces. So, I would need to know what quantities specifically does Abaqus provide and what quantities does it expect back from the model. From there, we can try and map those to what ExaCMech expects.

Also, the initialization of models is another thing to be mindful of as best practice from the ExaCMech side of things would be to only create the object once and then re-use it for the lifetime of things rather than recreating the object every material point call. We would also need to instantiate the state / history variables during the very first time step / first call per quadrature point as well to make sure all of the lattice orientations and history variables are instantiated correctly. Hopefully, Abaqus has a way for that to be done for UMATs. If not additional logic will be required to guard against such things.

I will also mention that I'm not opposed to adding these wrappers to ExaCMech as we can have them only built when an optional CMake option is defined, but I probably wouldn't be the one to create them due to time constraints. However, I would be able to do of mentoring type thing to help guide someone into creating that initial prototype / implementation of a wrapper.

vigneshwaran-radhakrishnan commented 4 days ago

@rcarson3 Thanks. At each step, UMAT provides the following data for every integration point:

  1. Initial data: All user-specified parameters
  2. Previous step data: Cauchy stress, log strain, deformation gradient, and state variables (including all user-defined slip system orientations, hardening parameters, etc..)
  3. Current step data: Deformation (strain) rate (D) and rotation rate

UMAT expects to return the updated Cauchy stress and the Jacobian matrix (del Stress/del Strain), and that's it.

I can put effort and time on to it. Guidance would be great. What is a good starting point? Start with calling pyecmech (or the jax version) fior every integration point from umat or something else?

rcarson3 commented 4 days ago

You probably don't want to use the python bindings as that will be slower. Plus, the python stuff is more for prototyping new models or running material point simulation type things to get a better understanding of a material models behavior.

You can take a look at the miniapp we have for the library: https://github.com/LLNL/ExaCMech/blob/develop/miniapp/ and more specifically: https://github.com/LLNL/ExaCMech/blob/develop/miniapp/orientation_evolution.cxx and it's subsequent calls to see how an application code might call the library and what data it will provide the code. So, I'd say understanding the differences in the calls there to what UMAT expects would be a good first step.

I will admit the miniapp isn't the cleanest code right now, but it was meant to provide people with a model for how to use the library at a higher level.

Within the miniapp itself, you'll note there is an init_data call which setups all the state variables and this is only called once. In the orientation_evolution file, you'll notice we only create the material model once as we don't need to recreate our c++ classes constantly as all of their member variables stay the same through out the life of the program. For each time step, we have 3 steps we go through, a setup phase called setup_data found in the setup_kernels.cxx file, a batched material model invocation called from material_kernels.cxx and then a post-processing stage called retrieve_data in the retrieve_kernels.cxx file.

If I recall UMATs are per quadrature point calls right and are not batch quadrature point calls right? If so we'd probably just need to tell the underlying ECMech calls that it only needs to do 1 point for the calculation which is sub-optimal from a performance point of view but should be fine for now.

If Abaqus provides which newton iteration we're on for the nonlinear solve and we know the time step is 0 then that should be enough info to properly set-up the data only once ever. Since, you can put guards into the code to check for such things. I'm guessing the instantiation of the crystal orientation as it varies from element to element is an already understood solution there so I'm guessing that's something y'all already know how to handle. Outside of that, yeah I think what I mentioned up above might be the best 1st steps and we can go from there.