Lithiation-delithiation#
Code for Cahn-Hilliard phase separation with mechanical coupling
This is a two-dimensional axisymmetric simulation
Slow ramp lithiation-delithiation
Degrees of freedom#
mechanical dispalcement: u
chemical potential: mu
concentration: c
Units#
Length: um
Mass: kg
Time: s
Amount of substance: pmol
Temperature: K
Mass density: kg/um^3
Force: uN
Stress: MPa
Energy: pJ
Species concentration: pmol/um^3
Chemical potential: pJ/pmol
Molar volume: um^3/pmol
Species diffusivity: um^2/s
Boltzmann Constant: 1.38E-11 pJ/K
Gas constant: 8.3145 pJ/(pmol K)
Software:#
Dolfinx v0.8.0
In the collection “Example Codes for Coupled Theories in Solid Mechanics,”
By Eric M. Stewart, Shawn A. Chester, and Lallit Anand.
Import modules#
# Import FEnicSx/dolfinx
import dolfinx
# For numerical arrays
import numpy as np
# For MPI-based parallelization
from mpi4py import MPI
comm = MPI.COMM_WORLD
rank = comm.Get_rank()
# PETSc solvers
from petsc4py import PETSc
# specific functions from dolfinx modules
from dolfinx import fem, mesh, io, plot, log
from dolfinx.fem import (Constant, dirichletbc, Function, functionspace, Expression )
from dolfinx.fem.petsc import NonlinearProblem
from dolfinx.nls.petsc import NewtonSolver
from dolfinx.io import VTXWriter, XDMFFile
# specific functions from ufl modules
import ufl
from ufl import (TestFunctions, TrialFunction, Identity, grad, det, div, dev, inv, tr, sqrt, conditional ,\
lt, gt, dx, inner, derivative, dot, ln, split, exp, eq, cos, sin, acos, ge, le, outer, tanh,\
cosh, atan, atan2)
# basix finite elements (necessary for dolfinx v0.8.0)
import basix
from basix.ufl import element, mixed_element
# Matplotlib for plotting
import matplotlib.pyplot as plt
plt.close('all')
# For timing the code
from datetime import datetime
# Set level of detail for log messages (integer)
# Guide:
# CRITICAL = 50, // errors that may lead to data corruption
# ERROR = 40, // things that HAVE gone wrong
# WARNING = 30, // things that MAY go wrong later
# INFO = 20, // information of general interest (includes solver info)
# PROGRESS = 16, // what's happening (broadly)
# TRACE = 13, // what's happening (in detail)
# DBG = 10 // sundry
#
log.set_log_level(log.LogLevel.WARNING)
Define geometry#
# Ellipse major axis
a0 = 0.5 # um
# Ellipse minor axis
b0 = 0.3 #um
# Read in the 2D mesh and cell tags
with XDMFFile(MPI.COMM_WORLD,"meshes/axi_ellipse.xdmf",'r') as infile:
domain = infile.read_mesh(name="Grid",xpath="/Xdmf/Domain")
cell_tags = infile.read_meshtags(domain,name="Grid")
domain.topology.create_connectivity(domain.topology.dim, domain.topology.dim-1)
# Also read in 1D facets for applying BCs
with XDMFFile(MPI.COMM_WORLD,"meshes/facet_axi_ellipse.xdmf",'r') as infile:
facet_tags = infile.read_meshtags(domain,name="Grid")
# A single point at the center of the inclusion for "grounding" the displacement of the inclusion
def ground(x):
return np.logical_and(np.isclose(x[0], 0), np.isclose(x[1], 0))
x = ufl.SpatialCoordinate(domain)
Print out the unique cell number indices
top_imap = domain.topology.index_map(2) # index map of 2D entities in domain
values = np.zeros(top_imap.size_global) # an array of zeros of the same size as number of 2D entities
values[cell_tags.indices]=cell_tags.values # populating the array with facet tag index numbers
print(np.unique(cell_tags.values)) # printing the unique indices
# Gmsh numbering
# Physical Surface("particle", 8)
[8]
Print out the unique edge index numbers
top_imap = domain.topology.index_map(1) # index map of 1D entities in domain
values = np.zeros(top_imap.size_global) # an array of zeros of the same size as number of 2D entities
values[facet_tags.indices]=facet_tags.values # populating the array with facet tag index numbers
print(np.unique(facet_tags.values)) # printing the unique indices
# Curve labels from Gmsh
#
# //+
# Physical Curve("left", 5)
# //+
# Physical Curve("bottom", 6)
# //+
# Physical Curve("outer", 7)
[5 6 7]
Visualize the reference configuration
import pyvista
pyvista.set_jupyter_backend('html')
from dolfinx.plot import vtk_mesh
pyvista.start_xvfb()
# initialize a plotter
plotter = pyvista.Plotter()
# Add the mesh.
topology, cell_types, geometry = plot.vtk_mesh(domain, domain.topology.dim)
grid = pyvista.UnstructuredGrid(topology, cell_types, geometry)
plotter.add_mesh(grid, show_edges=True) #, opacity=0.25)
plotter.view_xy()
#labels = dict(xlabel='X', ylabel='Y',zlabel='Z')
labels = dict(xlabel='X', ylabel='Y')
plotter.add_axes(**labels)
plotter.screenshot("results/ellipse_mesh.png")
from IPython.display import Image
Image(filename='results/ellipse_mesh.png')
# # Use the following commands for a zoom-able view
# if not pyvista.OFF_SCREEN:
# plotter.show()
# else:
# plotter.screenshot("square_mesh.png")

Define boundary and volume integration measure#
# Define the boundary integration measure "ds" using the facet tags,
# also specify the number of surface quadrature points.
ds = ufl.Measure('ds', domain=domain, subdomain_data=facet_tags, metadata={'quadrature_degree':4})
# Define the volume integration measure "dx"
# also specify the number of volume quadrature points.
dx = ufl.Measure('dx', domain=domain, metadata={'quadrature_degree': 4})
# Create facet to cell connectivity required to determine boundary facets.
domain.topology.create_connectivity(domain.topology.dim, domain.topology.dim)
domain.topology.create_connectivity(domain.topology.dim, domain.topology.dim-1)
domain.topology.create_connectivity(domain.topology.dim-1, domain.topology.dim)
# # Define facet normal
n2D = ufl.FacetNormal(domain)
n = ufl.as_vector([n2D[0], n2D[1], 0.0]) # define n as a 3D vector for later use
Material parmeters#
# Material parameters after Di Leo et al. (2014)
factor = Constant(domain,PETSc.ScalarType(1.0))
#
Gshear = Constant(domain,PETSc.ScalarType(49.8e3*factor)) # Shear modulus, MPa
Kbulk = Constant(domain,PETSc.ScalarType(83e3*factor)) # Bulk modulus, MPa
#
Omega = Constant(domain,PETSc.ScalarType(4.05)) # Molar volume, um^3/pmol
D = Constant(domain,PETSc.ScalarType(1e-2)) # Diffusivity, um^2/s
chi = Constant(domain,PETSc.ScalarType(3)) # Phase parameter, (-)
cMax = Constant(domain,PETSc.ScalarType(2.29e-2)) # Saturation concentration, pmol/um^3
lam = Constant(domain,PETSc.ScalarType(5.5749e-1)) # Interface parameter, (pJ/pmol) um^2
#
theta0 = Constant(domain,PETSc.ScalarType(298)) # Reference temperature, K
R_gas = Constant(domain,PETSc.ScalarType(8.3145)) # Gas constant, pJ/(pmol K)
RT = Constant(domain,PETSc.ScalarType(R_gas*theta0))
Function spaces#
# Define function space, both vectorial and scalar
U2 = element("Lagrange", domain.basix_cell(), 2, shape=(2,)) # For displacement
P1 = element("Lagrange", domain.basix_cell(), 1) # Forchemical potential and species concentration
#
TH = mixed_element([U2, P1, P1]) # Mixed element
ME = functionspace(domain, TH) # Total space for all DOFs
# Define actual functions with the required DOFs
w = Function(ME)
u, mu, c = split(w) # chemical potential mu and concentration c
# A copy of functions to store values in the previous step for time-stepping
w_old = Function(ME)
u_old, mu_old, c_old = split(w_old)
# Define test functions
u_test, mu_test, c_test = TestFunctions(ME)
# Define trial functions needed for automatic differentiation
dw = TrialFunction(ME)
Initial conditions#
# Assign initial normalized chemical potential mu_i to the domain
w.sub(1).interpolate(lambda x: np.full((x.shape[1],), mu_i))
w_old.sub(1).interpolate(lambda x: np.full((x.shape[1],), mu_i))
# Assign initial normalized concentration c_i to the domain
c_i = Constant(domain, PETSc.ScalarType(cBar_i*cMax*Omega))
#
w.sub(2).interpolate(lambda x: np.full((x.shape[1],), c_i))
w_old.sub(2).interpolate(lambda x: np.full((x.shape[1],), c_i))
Subroutines for kinematics and constitutive equations#
# Special gradient operators for axisymmetric functions
#
#---------------------------------------------------
# Gradient of vector field u
#---------------------------------------------------
def axi_grad_vector(u):
grad_u = grad(u)
axi_grad_33_exp = conditional(eq(x[0], 0), 0.0, u[0]/x[0])
axi_grad_u = ufl.as_tensor([[grad_u[0,0], grad_u[0,1], 0],
[grad_u[1,0], grad_u[1,1], 0],
[0, 0, axi_grad_33_exp]])
return axi_grad_u
#---------------------------------------------------
# Gradient of scalar field y
# (just need an extra zero for dimensions to work out)
#---------------------------------------------------
def axi_grad_scalar(y):
grad_y = grad(y)
axi_grad_y = ufl.as_vector([grad_y[0], grad_y[1], 0.])
return axi_grad_y
#---------------------------------------------------
# Axisymmetric deformation gradient
#---------------------------------------------------
def F_axi_calc(u):
dim = len(u) # dimension of problem (2)
Id = Identity(dim) # 2D Identity tensor
F = Id + grad(u) # 2D Deformation gradient
F33_exp = 1.0 + u[0]/x[0] # axisymmetric F33, R/R0
F33 = conditional(eq(x[0], 0.0), 1.0, F33_exp) # avoid divide by zero at r=0
F_axi = ufl.as_tensor([[F[0,0], F[0,1], 0],
[F[1,0], F[1,1], 0],
[0, 0, F33]]) # Full axisymmetric F
return F_axi
# Elastic deformation gradient Fe
def Fe_calc(u,c):
F = F_axi_calc(u) # = F
J = det(F) # = J
#
Js = 1.0 + c # what about (c-c_i)?
Fs = Js**(1/3)*Identity(3)
#
Fe = F*inv(Fs)
return Fe
# The elastic second Piola stress
def Te_calc(u, c):
Id = Identity(3)
#
Fe = Fe_calc(u, c)
Je = det(Fe)
Ce = Fe.T*Fe
#
Cebar = Je**(-2/3)*Ce
#
Te = Je**(-2/3)*Gshear*(Id - (1/3)*tr(Cebar)*inv(Cebar))\
+ Kbulk*Je*(Je-1)*inv(Cebar)
#
return Te
# The elastic Mandel stress
def Me_calc(u, c):
Fe = Fe_calc(u, c)
Je = det(Fe)
Ce = Fe.T*Fe
#
Te = Te_calc(u, c)
#
Me = Ce*Te
#
return Me
# The first Piola stress
def Piola_calc(u,c):
#
F = F_axi_calc(u)
J = det(F)
#
Fe = Fe_calc(u, c)
Je = det(Fe)
#
Te = Te_calc(u,c)
#
T = Je**(-1)*Fe*Te*inv(Fe)
#
Piola = J*T*inv(F.T)/Gshear
return Piola
#------------------------------------------------------------------------------
# Species flux
def Flux_calc(u, mu, c):
F = F_axi_calc(u)
#
Cinv = inv(F.T*F)
#
cBar = c/(Omega*cMax) # normalized concentration
#
Mob = (D*c)/(Omega*RT)*(1-cBar)*Cinv
#
Jmat = - RT* Mob * axi_grad_scalar(mu)
#
return Jmat
# Calculate the f^c term
def fc_calc(u, c):
#
cBar = c/(Omega*cMax) # normalized concentration
#
Me = Me_calc(u,c)
#
fc = RT*(ln(cBar/(1-cBar)) + chi*(1-2*cBar) ) - Omega*((1/3)*tr(Me))
#
return fc
# Calculate principal Cauchy stresses for visualization only
def tensor_eigs(T):
# invariants of T
I1 = tr(T)
I2 = (1/2)*(tr(T)**2 - tr(T*T))
I3 = det(T)
# Intermediate quantities b, c, d
b = -I1
c = I2
d = -I3
# intermediate quantities E, F, G
E = (3*c - b*b)/3
F = (2*(b**3) - 9*b*c + 27*d)/27
G = (F**2)/4 + (E**3)/27
# Intermediate quantities H, I, J, K, L
H = sqrt(-(E**3)/27)
I = H**(1/3)
J = acos(-F/(2*H))
K = cos(J/3)
L = sqrt(3)*sin(J/3)
# Finally, the (not necessarily ordered) eigenvalues
t1 = 2*I*K - b/3
t2 = -I*(K+L) - b/3
t3 = -I*(K-L) - b/3
# Order the eigenvalues using conditionals
#
T1_temp = conditional(lt(t1, t3), t3, t1 ) # returns the larger of t1 and t3.
T1 = conditional(lt(T1_temp, t2), t2, T1_temp ) # returns the larger of T1_temp and t2.
#
T3_temp = conditional(gt(t3, t1), t1, t3 ) # returns the smaller of t1 and t3.
T3 = conditional(gt(T3_temp, t2), t2, T1_temp ) # returns the smaller of T3_temp and t2.
#
# use the trace to report the middle eigenvalue.
T2 = I1 - T1 - T3
return T1, T2, T3
Evaluate kinematics and constitutive relations#
# Kinematics
F = F_axi_calc(u)
J = det(F)
# Calculate the normalized concentration cBar
cBar = c/(Omega*cMax) # normalized concentration
# Calculate the Piola stress
Piola = Piola_calc(u,c)
# Calculate the Species flux
Jmat = Flux_calc(u, mu, c)
# Calculate the f^c term
fc = fc_calc(u, c)
Weak forms#
# Residuals:
# Res_0: Equation of motion (test fxn: u)
# Res_1: Balance of mass (test fxn: mu)
# Res_2: chemical potential (test fxn: c)
#
mu_ext_cons = Constant(domain,PETSc.ScalarType(muRamp(0)))
# Define boundary flux :
j_flux = mu - mu_ext_cons
# # # We need to work on this:
# # # Calculate the spatial species flux:
# j_spat = mu - mu_cons
# #
# Fcof = J*inv(F.T) # Cofactor of F
# nvec = dot(Fcof,n) # Recall that we had defined: n = ufl.as_vector([n2D[0], n2D[1], 0.0])
# da_mat = sqrt(inner(nvec, nvec)) # Areal jacobian
# # Calculate the configuration-dependent referential spcies in-flux (note the minus sign)
# jBar = j_flux * da_mat
# The weak form for the equation of motion
Res_0 = inner(Piola, axi_grad_vector(u_test))*x[0]*dx
# The weak form for the mass balance of mobile species
Res_1 = dot((c - c_old)/dk, mu_test)*x[0]*dx \
- Omega*dot(Jmat , axi_grad_scalar(mu_test) )*x[0]*dx \
+ Omega*dot(j_flux, mu_test)*x[0]*ds(7)
# The weak form for the mass balance of mobile species
# Res_1 = dot((c - c_old)/dk, mu_test)*x[0]*dx \
# - Omega*dot(Jmat , axi_grad_scalar(mu_test) )*x[0]*dx \
# + Omega*dot(jBar, mu_test)*x[0]*ds(7)
# The weak form for the concentration
Res_2 = dot(mu - fc/RT, c_test)*x[0]*dx \
- dot( (lam/RT)*axi_grad_scalar(cBar), axi_grad_scalar(c_test))*x[0]*dx
# Total weak form
Res = Res_0 + Res_1 + Res_2
# Automatic differentiation tangent:
a = derivative(Res, w, dw)
Set-up output files#
# Set up projection problem for fixing visualization issues
# of fields in the axisymmetric simulation
#
def setup_projection(u, V):
trial = ufl.TrialFunction(V)
test = ufl.TestFunction(V)
a = ufl.inner(trial, test)*x[0]*dx
L = ufl.inner(u, test)*x[0]*dx
projection_problem = dolfinx.fem.petsc.LinearProblem(a, L, [], \
petsc_options={"ksp_type": "cg", "ksp_rtol": 1e-16, "ksp_atol": 1e-16, "ksp_max_it": 1000})
return projection_problem
# results file name
results_name = "axi_ellipse_lithiation_delithiation_mod_factor=1.0"
# Function space for projection of results
U1 = element("DG", domain.basix_cell(), 1, shape=(2,)) # For 2d vector
P0 = element("DG", domain.basix_cell(), 1) # For scalar
T1 = element("DG", domain.basix_cell(), 1, shape=(3,3)) # For stress tensor
V1 = fem.functionspace(domain, P0) # Scalar function space
V2 = fem.functionspace(domain, U1) # Vector function space
V3 = fem.functionspace(domain, T1) # Tensor function space
# basic fields to write to output file
u_vis = Function(V2)
u_vis.name = "disp"
mu_vis = Function(V1)
mu_vis.name = "mu"
c_vis = Function(V1)
c_vis.name = "c"
# Project cBar for visualization
#
cBar_projection = setup_projection(cBar, V1)
cBar_vis = cBar_projection.solve()
cBar_vis.name = "cBar"
# Project the volumetric Jacobian J for visualization
#
J_projection = setup_projection(J, V1)
J_vis = J_projection.solve()
J_vis.name = "J"
# Project the Piola stress tensor for visualization
#
Piola_projection = setup_projection(Piola, V3)
Piola_temp = Piola_projection.solve()
T = Gshear*Piola_temp*F.T/J
T0 = T - (1/3)*tr(T)*Identity(3)
Mises = sqrt((3/2)*inner(T0, T0))
Mises_projection = setup_projection(Mises, V1)
Mises_vis = Mises_projection.solve()
Mises_vis.name = "Mises"
sig1, sig2, sig3 = tensor_eigs(T)
#
sig1_projection = setup_projection(sig1, V1)
sig1_vis = sig1_projection.solve()
sig1_vis.name = "sig1"
#
sig2_projection = setup_projection(sig2, V1)
sig2_vis = sig2_projection.solve()
sig2_vis.name = "sig2"
#
sig3_projection = setup_projection(sig3, V1)
sig3_vis = sig3_projection.solve()
sig3_vis.name = "sig3"
# #
# sig1_vis = Function(V1)
# sig1_vis.name = "sig1"
# sig1_expr = Expression(sig1, V1.element.interpolation_points())
# #
# sig2_vis = Function(V1)
# sig2_vis.name = "sig2"
# sig2_expr = Expression(sig2, V1.element.interpolation_points())
P11 = Function(V1)
P11.name = "P11"
P11_expr = Expression(Gshear*Piola_temp[0,0],V1.element.interpolation_points())
#
P22 = Function(V1)
P22.name = "P22"
P22_expr = Expression(Gshear*Piola_temp[1,1],V1.element.interpolation_points())
#
P33 = Function(V1)
P33.name = "P33"
P33_expr = Expression(Gshear*Piola_temp[2,2],V1.element.interpolation_points())
# set up the output VTX files.
file_results = VTXWriter(
MPI.COMM_WORLD,
"results/" + results_name + ".bp",
[ # put the functions here you wish to write to output
u_vis, mu_vis, c_vis, cBar_vis, J_vis,
P11, P22, P33, sig1_vis, sig2_vis, sig3_vis, Mises_vis,
],
engine="BP4",
)
def writeResults(t):
# Re-project some fields. This is necessary here to remove visual artifacts which arise
# due to the axisymmetric formulation as r -> 0
#
Piola_temp = Piola_projection.solve()
Mises_vis = Mises_projection.solve()
J_vis = J_projection.solve()
cBar_vis = cBar_projection.solve()
sig1_vis = sig1_projection.solve()
sig2_vis = sig2_projection.solve()
sig3_vis = sig3_projection.solve()
# Output field interpolation
u_vis.interpolate(w.sub(0))
mu_vis.interpolate(w.sub(1))
c_vis.interpolate(w.sub(2))
P11.interpolate(P11_expr)
P22.interpolate(P22_expr)
P33.interpolate(P33_expr)
# #
# sig1_vis.interpolate(sig1_expr)
# sig2_vis.interpolate(sig2_expr)
# Write output fields
file_results.write(t)
Infrastructure for pulling out time history data (displacement, force, etc.)#
# Identify point for reporting the chemical potential.
# The coordinates belwo are for the outer point of the minor axis of the ellipse
pointForDisp = np.array([b0,0.0,0.0])
bb_tree = dolfinx.geometry.bb_tree(domain,domain.topology.dim)
cell_candidates = dolfinx.geometry.compute_collisions_points(bb_tree, pointForDisp)
colliding_cells = dolfinx.geometry.compute_colliding_cells(domain, cell_candidates, pointForDisp).array
# Form for evaluating the state of charge SOC at each step
# ( SOC = Tot_chg / Tot_vol )
#
Tot_chg = fem.form(2.0*np.pi* cBar *x[0]*dx)
#
Tot_vol = fem.form(2.0*np.pi*x[0]*dx)
# The (referential) volume is constant, so just evaluate it once before running the simulation:
#
Tot_vol_val = domain.comm.gather(fem.assemble_scalar(Tot_vol))[0]
Analysis Step#
# Give the step a descriptive name
step = "Lithiate"
Boundary conditions#
# Constant for applied chemical potential
mu_cons = Constant(domain,PETSc.ScalarType(muRamp(0)))
# Recall gmsh curve names and numbers
# Physical Curve("left", 5)
# //+
# Physical Curve("bottom", 6)
# //+
# Physical Curve("outer", 7)
# Find the specific DOFs which will be constrained
xBot_u1_dofs = fem.locate_dofs_topological(ME.sub(0).sub(0), facet_tags.dim, facet_tags.find(5))
yBot_u2_dofs = fem.locate_dofs_topological(ME.sub(0).sub(1), facet_tags.dim, facet_tags.find(6))
#
#outer_mu_dofs = fem.locate_dofs_topological(ME.sub(1), facet_tags.dim, facet_tags.find(7))
# Build Dirichlet BCs
bcs_1 = dirichletbc(0.0, xBot_u1_dofs, ME.sub(0).sub(0)) # u1 fix - left
bcs_2 = dirichletbc(0.0, yBot_u2_dofs, ME.sub(0).sub(1)) # u2 fix - bottom
#
bcs = [bcs_1,bcs_2]
Define the nonlinear variational problem#
# Set up nonlinear problem
problem = NonlinearProblem(Res, w, bcs, a)
# The global newton solver and params
solver = NewtonSolver(MPI.COMM_WORLD, problem)
solver.convergence_criterion = "incremental"
solver.rtol = 1e-8
solver.atol = 1e-8
solver.max_it = 50
solver.report = True
solver.error_on_nonconvergence = False
# The Krylov solver parameters.
ksp = solver.krylov_solver
opts = PETSc.Options()
option_prefix = ksp.getOptionsPrefix()
opts[f"{option_prefix}ksp_type"] = "preonly"
opts[f"{option_prefix}pc_type"] = "lu" # do not use 'gamg' pre-conditioner
opts[f"{option_prefix}pc_factor_mat_solver_type"] = "mumps"
opts[f"{option_prefix}ksp_max_it"] = 30
ksp.setFromOptions()
Initialize arrays for storing output history#
#Arrays for storing output history
totSteps = 100000
timeHist0 = np.zeros(shape=[totSteps])
timeHist1 = np.zeros(shape=[totSteps])
timeHist2 = np.zeros(shape=[totSteps])
# # timeHist3 = np.zeros(shape=[totSteps])
# # #
timeHist1[0] = mu_i # Initialize the chemical potential
# Initialize a counter for reporting data
ii = 0
# Write initial state to file
writeResults(t=0.0)
Start calculation loop#
print("------------------------------------")
print("Simulation start")
print("------------------------------------")
# Store start time
startTime = datetime.now()
# Time-stepping solution procedure loop
while (round(t, 9) <= Ttot):
# increment time
t += dt
# update time variables in time-dependent BCs
mu_ext_cons.value = float(muRamp(t))
# Solve the problem
(iter, converged) = solver.solve(w)
# Now we start the adaptive time-stepping and output storage procedure.
#
# First, we check if the newton solver actually converged.
if converged:
# Collect results from MPI ghost processes
w.x.scatter_forward()
# Write output to file
writeResults(t)
# If the solver converged, we print the status of the solver,
# perform adaptive time-stepping updates, output results, and
# update degrees of freedom for the next step, w_old <- w.
# print progress of calculation periodically
if ii%1 == 0:
now = datetime.now()
current_time = now.strftime("%H:%M:%S")
print("Step: {} | Increment: {} | Iterations: {}".format(step, ii, iter))
print("Simulation Time: {} s | dt: {} s".format(round(t,2), round(dt, 3)))
print()
# Iteration-based adaptive time-stepping
#
# If the newton solver takes 2 or less iterations,
# increase the time step by a factor of 1.5:
if ((iter<=4) and (dt<25)):
dt = 1.5*dt
dk.value = dt
# If the newton solver takes 5 or more iterations,
# decrease the time step by a factor of 2:
elif iter>=8:
dt = dt/2
dk.value = dt
# otherwise the newton solver took 3-4 iterations,
# in which case leave the time step alone.
# update the DOFs for the next step.
w_old.x.array[:] = w.x.array
# Increment counter
ii += 1
# Store time history variables
timeHist0[ii] = t # current time
timeHist1[ii] = w.sub(1).eval([b0,0.0,0.0],colliding_cells[0])[0] # time history of chemical potetial at aouter boundary
timeHist2[ii] = domain.comm.gather(fem.assemble_scalar(Tot_chg))[0] / Tot_vol_val # time history of State of charge
# If solver doesn't converge we have to back up in time,
# cut the size of the time step, and try solving again.
else: # not(converged)
# first, we back up in time
# ( to un-do the current time step )
t = t - float(dk)
# Then, we cut back on the time step we're attempting.
# (by a factor of 2)
dt = dt/2
dk.value = dt
# Re-set the DOFs to their value before the failed step.
w.x.array[:] = w_old.x.array
# Set up nonlinear problem
problem = NonlinearProblem(Res, w, bcs, a)
# The global newton solver and params
solver = NewtonSolver(MPI.COMM_WORLD, problem)
solver.convergence_criterion = "incremental"
solver.rtol = 1e-8
solver.atol = 1e-8
solver.max_it = 50
solver.report = True
solver.error_on_nonconvergence = False
# close the output file.
file_results.close()
# End analysis
print("-----------------------------------------")
print("End computation")
# Report elapsed real time for the analysis
endTime = datetime.now()
elapseTime = endTime - startTime
print("------------------------------------------")
print("Elapsed real time: {}".format(elapseTime))
print("------------------------------------------")
------------------------------------
Simulation start
------------------------------------
Step: Lithiate | Increment: 0 | Iterations: 3
Simulation Time: 1.0 s | dt: 1.0 s
Step: Lithiate | Increment: 1 | Iterations: 3
Simulation Time: 2.5 s | dt: 1.5 s
Step: Lithiate | Increment: 2 | Iterations: 3
Simulation Time: 4.75 s | dt: 2.25 s
Step: Lithiate | Increment: 3 | Iterations: 3
Simulation Time: 8.12 s | dt: 3.375 s
Step: Lithiate | Increment: 4 | Iterations: 4
Simulation Time: 13.19 s | dt: 5.062 s
Step: Lithiate | Increment: 5 | Iterations: 4
Simulation Time: 20.78 s | dt: 7.594 s
Step: Lithiate | Increment: 6 | Iterations: 4
Simulation Time: 32.17 s | dt: 11.391 s
Step: Lithiate | Increment: 7 | Iterations: 4
Simulation Time: 49.26 s | dt: 17.086 s
Step: Lithiate | Increment: 8 | Iterations: 4
Simulation Time: 74.89 s | dt: 25.629 s
Step: Lithiate | Increment: 9 | Iterations: 4
Simulation Time: 100.52 s | dt: 25.629 s
Step: Lithiate | Increment: 10 | Iterations: 4
Simulation Time: 126.14 s | dt: 25.629 s
Step: Lithiate | Increment: 11 | Iterations: 4
Simulation Time: 151.77 s | dt: 25.629 s
Step: Lithiate | Increment: 12 | Iterations: 4
Simulation Time: 177.4 s | dt: 25.629 s
Step: Lithiate | Increment: 13 | Iterations: 4
Simulation Time: 203.03 s | dt: 25.629 s
Step: Lithiate | Increment: 14 | Iterations: 4
Simulation Time: 228.66 s | dt: 25.629 s
Step: Lithiate | Increment: 15 | Iterations: 4
Simulation Time: 254.29 s | dt: 25.629 s
Step: Lithiate | Increment: 16 | Iterations: 4
Simulation Time: 279.92 s | dt: 25.629 s
Step: Lithiate | Increment: 17 | Iterations: 4
Simulation Time: 305.55 s | dt: 25.629 s
Step: Lithiate | Increment: 18 | Iterations: 4
Simulation Time: 331.18 s | dt: 25.629 s
Step: Lithiate | Increment: 19 | Iterations: 4
Simulation Time: 356.8 s | dt: 25.629 s
Step: Lithiate | Increment: 20 | Iterations: 4
Simulation Time: 382.43 s | dt: 25.629 s
Step: Lithiate | Increment: 21 | Iterations: 4
Simulation Time: 408.06 s | dt: 25.629 s
Step: Lithiate | Increment: 22 | Iterations: 4
Simulation Time: 433.69 s | dt: 25.629 s
Step: Lithiate | Increment: 23 | Iterations: 4
Simulation Time: 459.32 s | dt: 25.629 s
Step: Lithiate | Increment: 24 | Iterations: 4
Simulation Time: 484.95 s | dt: 25.629 s
Step: Lithiate | Increment: 25 | Iterations: 4
Simulation Time: 510.58 s | dt: 25.629 s
Step: Lithiate | Increment: 26 | Iterations: 4
Simulation Time: 536.21 s | dt: 25.629 s
Step: Lithiate | Increment: 27 | Iterations: 4
Simulation Time: 561.84 s | dt: 25.629 s
Step: Lithiate | Increment: 28 | Iterations: 4
Simulation Time: 587.46 s | dt: 25.629 s
Step: Lithiate | Increment: 29 | Iterations: 4
Simulation Time: 613.09 s | dt: 25.629 s
Step: Lithiate | Increment: 30 | Iterations: 4
Simulation Time: 638.72 s | dt: 25.629 s
Step: Lithiate | Increment: 31 | Iterations: 4
Simulation Time: 664.35 s | dt: 25.629 s
Step: Lithiate | Increment: 32 | Iterations: 4
Simulation Time: 689.98 s | dt: 25.629 s
Step: Lithiate | Increment: 33 | Iterations: 4
Simulation Time: 715.61 s | dt: 25.629 s
Step: Lithiate | Increment: 34 | Iterations: 4
Simulation Time: 741.24 s | dt: 25.629 s
Step: Lithiate | Increment: 35 | Iterations: 4
Simulation Time: 766.87 s | dt: 25.629 s
Step: Lithiate | Increment: 36 | Iterations: 4
Simulation Time: 792.5 s | dt: 25.629 s
Step: Lithiate | Increment: 37 | Iterations: 4
Simulation Time: 818.12 s | dt: 25.629 s
Step: Lithiate | Increment: 38 | Iterations: 4
Simulation Time: 843.75 s | dt: 25.629 s
Step: Lithiate | Increment: 39 | Iterations: 4
Simulation Time: 869.38 s | dt: 25.629 s
Step: Lithiate | Increment: 40 | Iterations: 4
Simulation Time: 895.01 s | dt: 25.629 s
Step: Lithiate | Increment: 41 | Iterations: 4
Simulation Time: 920.64 s | dt: 25.629 s
Step: Lithiate | Increment: 42 | Iterations: 4
Simulation Time: 946.27 s | dt: 25.629 s
Step: Lithiate | Increment: 43 | Iterations: 4
Simulation Time: 971.9 s | dt: 25.629 s
Step: Lithiate | Increment: 44 | Iterations: 5
Simulation Time: 997.53 s | dt: 25.629 s
Step: Lithiate | Increment: 45 | Iterations: 5
Simulation Time: 1023.16 s | dt: 25.629 s
Step: Lithiate | Increment: 46 | Iterations: 5
Simulation Time: 1048.79 s | dt: 25.629 s
2025-02-13 16:18:19.446 ( 17.137s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 47 | Iterations: 6
Simulation Time: 1061.6 s | dt: 12.814 s
2025-02-13 16:18:20.654 ( 18.345s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
2025-02-13 16:18:22.136 ( 19.827s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 48 | Iterations: 5
Simulation Time: 1064.8 s | dt: 3.204 s
2025-02-13 16:18:24.666 ( 22.357s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 49 | Iterations: 5
Simulation Time: 1066.41 s | dt: 1.602 s
Step: Lithiate | Increment: 50 | Iterations: 6
Simulation Time: 1068.01 s | dt: 1.602 s
2025-02-13 16:18:26.142 ( 23.833s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
2025-02-13 16:18:27.217 ( 24.908s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 51 | Iterations: 4
Simulation Time: 1068.41 s | dt: 0.4 s
2025-02-13 16:18:28.654 ( 26.345s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 52 | Iterations: 4
Simulation Time: 1068.71 s | dt: 0.3 s
2025-02-13 16:18:29.882 ( 27.573s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 53 | Iterations: 4
Simulation Time: 1068.93 s | dt: 0.225 s
2025-02-13 16:18:31.117 ( 28.808s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 54 | Iterations: 4
Simulation Time: 1069.1 s | dt: 0.169 s
2025-02-13 16:18:32.267 ( 29.958s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 55 | Iterations: 4
Simulation Time: 1069.23 s | dt: 0.127 s
Step: Lithiate | Increment: 56 | Iterations: 7
Simulation Time: 1069.42 s | dt: 0.19 s
2025-02-13 16:18:33.694 ( 31.385s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 57 | Iterations: 6
Simulation Time: 1069.51 s | dt: 0.095 s
Step: Lithiate | Increment: 58 | Iterations: 6
Simulation Time: 1069.61 s | dt: 0.095 s
Step: Lithiate | Increment: 59 | Iterations: 6
Simulation Time: 1069.7 s | dt: 0.095 s
Step: Lithiate | Increment: 60 | Iterations: 6
Simulation Time: 1069.8 s | dt: 0.095 s
Step: Lithiate | Increment: 61 | Iterations: 5
Simulation Time: 1069.89 s | dt: 0.095 s
Step: Lithiate | Increment: 62 | Iterations: 5
Simulation Time: 1069.99 s | dt: 0.095 s
Step: Lithiate | Increment: 63 | Iterations: 5
Simulation Time: 1070.08 s | dt: 0.095 s
Step: Lithiate | Increment: 64 | Iterations: 5
Simulation Time: 1070.18 s | dt: 0.095 s
Step: Lithiate | Increment: 65 | Iterations: 5
Simulation Time: 1070.27 s | dt: 0.095 s
Step: Lithiate | Increment: 66 | Iterations: 5
Simulation Time: 1070.37 s | dt: 0.095 s
Step: Lithiate | Increment: 67 | Iterations: 5
Simulation Time: 1070.46 s | dt: 0.095 s
Step: Lithiate | Increment: 68 | Iterations: 5
Simulation Time: 1070.56 s | dt: 0.095 s
Step: Lithiate | Increment: 69 | Iterations: 5
Simulation Time: 1070.65 s | dt: 0.095 s
Step: Lithiate | Increment: 70 | Iterations: 5
Simulation Time: 1070.75 s | dt: 0.095 s
Step: Lithiate | Increment: 71 | Iterations: 5
Simulation Time: 1070.84 s | dt: 0.095 s
Step: Lithiate | Increment: 72 | Iterations: 5
Simulation Time: 1070.94 s | dt: 0.095 s
Step: Lithiate | Increment: 73 | Iterations: 5
Simulation Time: 1071.03 s | dt: 0.095 s
Step: Lithiate | Increment: 74 | Iterations: 5
Simulation Time: 1071.13 s | dt: 0.095 s
Step: Lithiate | Increment: 75 | Iterations: 5
Simulation Time: 1071.22 s | dt: 0.095 s
Step: Lithiate | Increment: 76 | Iterations: 5
Simulation Time: 1071.32 s | dt: 0.095 s
Step: Lithiate | Increment: 77 | Iterations: 5
Simulation Time: 1071.41 s | dt: 0.095 s
Step: Lithiate | Increment: 78 | Iterations: 5
Simulation Time: 1071.51 s | dt: 0.095 s
Step: Lithiate | Increment: 79 | Iterations: 5
Simulation Time: 1071.6 s | dt: 0.095 s
Step: Lithiate | Increment: 80 | Iterations: 5
Simulation Time: 1071.7 s | dt: 0.095 s
Step: Lithiate | Increment: 81 | Iterations: 5
Simulation Time: 1071.79 s | dt: 0.095 s
Step: Lithiate | Increment: 82 | Iterations: 5
Simulation Time: 1071.89 s | dt: 0.095 s
Step: Lithiate | Increment: 83 | Iterations: 5
Simulation Time: 1071.98 s | dt: 0.095 s
Step: Lithiate | Increment: 84 | Iterations: 5
Simulation Time: 1072.08 s | dt: 0.095 s
Step: Lithiate | Increment: 85 | Iterations: 5
Simulation Time: 1072.17 s | dt: 0.095 s
Step: Lithiate | Increment: 86 | Iterations: 5
Simulation Time: 1072.27 s | dt: 0.095 s
Step: Lithiate | Increment: 87 | Iterations: 5
Simulation Time: 1072.36 s | dt: 0.095 s
Step: Lithiate | Increment: 88 | Iterations: 5
Simulation Time: 1072.46 s | dt: 0.095 s
Step: Lithiate | Increment: 89 | Iterations: 5
Simulation Time: 1072.55 s | dt: 0.095 s
Step: Lithiate | Increment: 90 | Iterations: 5
Simulation Time: 1072.65 s | dt: 0.095 s
Step: Lithiate | Increment: 91 | Iterations: 5
Simulation Time: 1072.74 s | dt: 0.095 s
Step: Lithiate | Increment: 92 | Iterations: 5
Simulation Time: 1072.84 s | dt: 0.095 s
Step: Lithiate | Increment: 93 | Iterations: 5
Simulation Time: 1072.93 s | dt: 0.095 s
Step: Lithiate | Increment: 94 | Iterations: 5
Simulation Time: 1073.03 s | dt: 0.095 s
Step: Lithiate | Increment: 95 | Iterations: 5
Simulation Time: 1073.12 s | dt: 0.095 s
Step: Lithiate | Increment: 96 | Iterations: 5
Simulation Time: 1073.22 s | dt: 0.095 s
Step: Lithiate | Increment: 97 | Iterations: 5
Simulation Time: 1073.31 s | dt: 0.095 s
Step: Lithiate | Increment: 98 | Iterations: 5
Simulation Time: 1073.41 s | dt: 0.095 s
Step: Lithiate | Increment: 99 | Iterations: 5
Simulation Time: 1073.5 s | dt: 0.095 s
Step: Lithiate | Increment: 100 | Iterations: 5
Simulation Time: 1073.6 s | dt: 0.095 s
Step: Lithiate | Increment: 101 | Iterations: 5
Simulation Time: 1073.69 s | dt: 0.095 s
Step: Lithiate | Increment: 102 | Iterations: 5
Simulation Time: 1073.79 s | dt: 0.095 s
Step: Lithiate | Increment: 103 | Iterations: 5
Simulation Time: 1073.88 s | dt: 0.095 s
Step: Lithiate | Increment: 104 | Iterations: 5
Simulation Time: 1073.98 s | dt: 0.095 s
Step: Lithiate | Increment: 105 | Iterations: 5
Simulation Time: 1074.08 s | dt: 0.095 s
Step: Lithiate | Increment: 106 | Iterations: 4
Simulation Time: 1074.17 s | dt: 0.095 s
Step: Lithiate | Increment: 107 | Iterations: 5
Simulation Time: 1074.31 s | dt: 0.143 s
Step: Lithiate | Increment: 108 | Iterations: 5
Simulation Time: 1074.46 s | dt: 0.143 s
Step: Lithiate | Increment: 109 | Iterations: 5
Simulation Time: 1074.6 s | dt: 0.143 s
Step: Lithiate | Increment: 110 | Iterations: 5
Simulation Time: 1074.74 s | dt: 0.143 s
Step: Lithiate | Increment: 111 | Iterations: 5
Simulation Time: 1074.88 s | dt: 0.143 s
Step: Lithiate | Increment: 112 | Iterations: 5
Simulation Time: 1075.03 s | dt: 0.143 s
Step: Lithiate | Increment: 113 | Iterations: 5
Simulation Time: 1075.17 s | dt: 0.143 s
Step: Lithiate | Increment: 114 | Iterations: 5
Simulation Time: 1075.31 s | dt: 0.143 s
Step: Lithiate | Increment: 115 | Iterations: 5
Simulation Time: 1075.45 s | dt: 0.143 s
Step: Lithiate | Increment: 116 | Iterations: 5
Simulation Time: 1075.6 s | dt: 0.143 s
Step: Lithiate | Increment: 117 | Iterations: 5
Simulation Time: 1075.74 s | dt: 0.143 s
Step: Lithiate | Increment: 118 | Iterations: 5
Simulation Time: 1075.88 s | dt: 0.143 s
Step: Lithiate | Increment: 119 | Iterations: 5
Simulation Time: 1076.02 s | dt: 0.143 s
Step: Lithiate | Increment: 120 | Iterations: 5
Simulation Time: 1076.17 s | dt: 0.143 s
Step: Lithiate | Increment: 121 | Iterations: 4
Simulation Time: 1076.31 s | dt: 0.143 s
Step: Lithiate | Increment: 122 | Iterations: 5
Simulation Time: 1076.52 s | dt: 0.214 s
Step: Lithiate | Increment: 123 | Iterations: 5
Simulation Time: 1076.74 s | dt: 0.214 s
Step: Lithiate | Increment: 124 | Iterations: 5
Simulation Time: 1076.95 s | dt: 0.214 s
Step: Lithiate | Increment: 125 | Iterations: 5
Simulation Time: 1077.16 s | dt: 0.214 s
Step: Lithiate | Increment: 126 | Iterations: 5
Simulation Time: 1077.38 s | dt: 0.214 s
Step: Lithiate | Increment: 127 | Iterations: 5
Simulation Time: 1077.59 s | dt: 0.214 s
Step: Lithiate | Increment: 128 | Iterations: 5
Simulation Time: 1077.8 s | dt: 0.214 s
Step: Lithiate | Increment: 129 | Iterations: 5
Simulation Time: 1078.02 s | dt: 0.214 s
Step: Lithiate | Increment: 130 | Iterations: 5
Simulation Time: 1078.23 s | dt: 0.214 s
Step: Lithiate | Increment: 131 | Iterations: 5
Simulation Time: 1078.45 s | dt: 0.214 s
Step: Lithiate | Increment: 132 | Iterations: 5
Simulation Time: 1078.66 s | dt: 0.214 s
Step: Lithiate | Increment: 133 | Iterations: 5
Simulation Time: 1078.87 s | dt: 0.214 s
Step: Lithiate | Increment: 134 | Iterations: 5
Simulation Time: 1079.09 s | dt: 0.214 s
Step: Lithiate | Increment: 135 | Iterations: 5
Simulation Time: 1079.3 s | dt: 0.214 s
Step: Lithiate | Increment: 136 | Iterations: 5
Simulation Time: 1079.52 s | dt: 0.214 s
Step: Lithiate | Increment: 137 | Iterations: 5
Simulation Time: 1079.73 s | dt: 0.214 s
Step: Lithiate | Increment: 138 | Iterations: 5
Simulation Time: 1079.94 s | dt: 0.214 s
Step: Lithiate | Increment: 139 | Iterations: 5
Simulation Time: 1080.16 s | dt: 0.214 s
Step: Lithiate | Increment: 140 | Iterations: 5
Simulation Time: 1080.37 s | dt: 0.214 s
Step: Lithiate | Increment: 141 | Iterations: 5
Simulation Time: 1080.58 s | dt: 0.214 s
Step: Lithiate | Increment: 142 | Iterations: 5
Simulation Time: 1080.8 s | dt: 0.214 s
Step: Lithiate | Increment: 143 | Iterations: 5
Simulation Time: 1081.01 s | dt: 0.214 s
Step: Lithiate | Increment: 144 | Iterations: 5
Simulation Time: 1081.23 s | dt: 0.214 s
Step: Lithiate | Increment: 145 | Iterations: 5
Simulation Time: 1081.44 s | dt: 0.214 s
Step: Lithiate | Increment: 146 | Iterations: 5
Simulation Time: 1081.65 s | dt: 0.214 s
Step: Lithiate | Increment: 147 | Iterations: 5
Simulation Time: 1081.87 s | dt: 0.214 s
Step: Lithiate | Increment: 148 | Iterations: 5
Simulation Time: 1082.08 s | dt: 0.214 s
Step: Lithiate | Increment: 149 | Iterations: 5
Simulation Time: 1082.3 s | dt: 0.214 s
Step: Lithiate | Increment: 150 | Iterations: 5
Simulation Time: 1082.51 s | dt: 0.214 s
Step: Lithiate | Increment: 151 | Iterations: 6
Simulation Time: 1082.72 s | dt: 0.214 s
Step: Lithiate | Increment: 152 | Iterations: 6
Simulation Time: 1082.94 s | dt: 0.214 s
Step: Lithiate | Increment: 153 | Iterations: 7
Simulation Time: 1083.15 s | dt: 0.214 s
Step: Lithiate | Increment: 154 | Iterations: 5
Simulation Time: 1083.36 s | dt: 0.214 s
Step: Lithiate | Increment: 155 | Iterations: 5
Simulation Time: 1083.58 s | dt: 0.214 s
Step: Lithiate | Increment: 156 | Iterations: 5
Simulation Time: 1083.79 s | dt: 0.214 s
Step: Lithiate | Increment: 157 | Iterations: 4
Simulation Time: 1084.01 s | dt: 0.214 s
Step: Lithiate | Increment: 158 | Iterations: 4
Simulation Time: 1084.33 s | dt: 0.321 s
Step: Lithiate | Increment: 159 | Iterations: 4
Simulation Time: 1084.81 s | dt: 0.481 s
Step: Lithiate | Increment: 160 | Iterations: 4
Simulation Time: 1085.53 s | dt: 0.722 s
Step: Lithiate | Increment: 161 | Iterations: 4
Simulation Time: 1086.61 s | dt: 1.082 s
Step: Lithiate | Increment: 162 | Iterations: 4
Simulation Time: 1088.24 s | dt: 1.624 s
Step: Lithiate | Increment: 163 | Iterations: 4
Simulation Time: 1090.67 s | dt: 2.435 s
Step: Lithiate | Increment: 164 | Iterations: 4
Simulation Time: 1094.32 s | dt: 3.653 s
Step: Lithiate | Increment: 165 | Iterations: 4
Simulation Time: 1099.8 s | dt: 5.48 s
Step: Lithiate | Increment: 166 | Iterations: 4
Simulation Time: 1108.02 s | dt: 8.22 s
Step: Lithiate | Increment: 167 | Iterations: 4
Simulation Time: 1120.35 s | dt: 12.33 s
Step: Lithiate | Increment: 168 | Iterations: 4
Simulation Time: 1138.85 s | dt: 18.495 s
Step: Lithiate | Increment: 169 | Iterations: 4
Simulation Time: 1166.59 s | dt: 27.742 s
Step: Lithiate | Increment: 170 | Iterations: 4
Simulation Time: 1194.33 s | dt: 27.742 s
Step: Lithiate | Increment: 171 | Iterations: 4
Simulation Time: 1222.07 s | dt: 27.742 s
Step: Lithiate | Increment: 172 | Iterations: 4
Simulation Time: 1249.81 s | dt: 27.742 s
Step: Lithiate | Increment: 173 | Iterations: 4
Simulation Time: 1277.56 s | dt: 27.742 s
Step: Lithiate | Increment: 174 | Iterations: 4
Simulation Time: 1305.3 s | dt: 27.742 s
Step: Lithiate | Increment: 175 | Iterations: 4
Simulation Time: 1333.04 s | dt: 27.742 s
Step: Lithiate | Increment: 176 | Iterations: 4
Simulation Time: 1360.78 s | dt: 27.742 s
Step: Lithiate | Increment: 177 | Iterations: 4
Simulation Time: 1388.52 s | dt: 27.742 s
Step: Lithiate | Increment: 178 | Iterations: 4
Simulation Time: 1416.27 s | dt: 27.742 s
Step: Lithiate | Increment: 179 | Iterations: 4
Simulation Time: 1444.01 s | dt: 27.742 s
Step: Lithiate | Increment: 180 | Iterations: 4
Simulation Time: 1471.75 s | dt: 27.742 s
Step: Lithiate | Increment: 181 | Iterations: 4
Simulation Time: 1499.49 s | dt: 27.742 s
Step: Lithiate | Increment: 182 | Iterations: 4
Simulation Time: 1527.23 s | dt: 27.742 s
Step: Lithiate | Increment: 183 | Iterations: 4
Simulation Time: 1554.97 s | dt: 27.742 s
Step: Lithiate | Increment: 184 | Iterations: 4
Simulation Time: 1582.72 s | dt: 27.742 s
Step: Lithiate | Increment: 185 | Iterations: 4
Simulation Time: 1610.46 s | dt: 27.742 s
Step: Lithiate | Increment: 186 | Iterations: 4
Simulation Time: 1638.2 s | dt: 27.742 s
Step: Lithiate | Increment: 187 | Iterations: 4
Simulation Time: 1665.94 s | dt: 27.742 s
Step: Lithiate | Increment: 188 | Iterations: 4
Simulation Time: 1693.68 s | dt: 27.742 s
Step: Lithiate | Increment: 189 | Iterations: 4
Simulation Time: 1721.42 s | dt: 27.742 s
Step: Lithiate | Increment: 190 | Iterations: 4
Simulation Time: 1749.17 s | dt: 27.742 s
Step: Lithiate | Increment: 191 | Iterations: 4
Simulation Time: 1776.91 s | dt: 27.742 s
Step: Lithiate | Increment: 192 | Iterations: 4
Simulation Time: 1804.65 s | dt: 27.742 s
Step: Lithiate | Increment: 193 | Iterations: 4
Simulation Time: 1832.39 s | dt: 27.742 s
Step: Lithiate | Increment: 194 | Iterations: 4
Simulation Time: 1860.13 s | dt: 27.742 s
Step: Lithiate | Increment: 195 | Iterations: 4
Simulation Time: 1887.88 s | dt: 27.742 s
Step: Lithiate | Increment: 196 | Iterations: 4
Simulation Time: 1915.62 s | dt: 27.742 s
Step: Lithiate | Increment: 197 | Iterations: 4
Simulation Time: 1943.36 s | dt: 27.742 s
Step: Lithiate | Increment: 198 | Iterations: 4
Simulation Time: 1971.1 s | dt: 27.742 s
Step: Lithiate | Increment: 199 | Iterations: 4
Simulation Time: 1998.84 s | dt: 27.742 s
Step: Lithiate | Increment: 200 | Iterations: 4
Simulation Time: 2026.58 s | dt: 27.742 s
Step: Lithiate | Increment: 201 | Iterations: 4
Simulation Time: 2054.33 s | dt: 27.742 s
Step: Lithiate | Increment: 202 | Iterations: 4
Simulation Time: 2082.07 s | dt: 27.742 s
Step: Lithiate | Increment: 203 | Iterations: 4
Simulation Time: 2109.81 s | dt: 27.742 s
Step: Lithiate | Increment: 204 | Iterations: 4
Simulation Time: 2137.55 s | dt: 27.742 s
Step: Lithiate | Increment: 205 | Iterations: 4
Simulation Time: 2165.29 s | dt: 27.742 s
Step: Lithiate | Increment: 206 | Iterations: 4
Simulation Time: 2193.03 s | dt: 27.742 s
Step: Lithiate | Increment: 207 | Iterations: 4
Simulation Time: 2220.78 s | dt: 27.742 s
Step: Lithiate | Increment: 208 | Iterations: 4
Simulation Time: 2248.52 s | dt: 27.742 s
Step: Lithiate | Increment: 209 | Iterations: 4
Simulation Time: 2276.26 s | dt: 27.742 s
Step: Lithiate | Increment: 210 | Iterations: 4
Simulation Time: 2304.0 s | dt: 27.742 s
Step: Lithiate | Increment: 211 | Iterations: 4
Simulation Time: 2331.74 s | dt: 27.742 s
Step: Lithiate | Increment: 212 | Iterations: 4
Simulation Time: 2359.49 s | dt: 27.742 s
Step: Lithiate | Increment: 213 | Iterations: 4
Simulation Time: 2387.23 s | dt: 27.742 s
Step: Lithiate | Increment: 214 | Iterations: 4
Simulation Time: 2414.97 s | dt: 27.742 s
Step: Lithiate | Increment: 215 | Iterations: 4
Simulation Time: 2442.71 s | dt: 27.742 s
Step: Lithiate | Increment: 216 | Iterations: 4
Simulation Time: 2470.45 s | dt: 27.742 s
Step: Lithiate | Increment: 217 | Iterations: 4
Simulation Time: 2498.19 s | dt: 27.742 s
Step: Lithiate | Increment: 218 | Iterations: 4
Simulation Time: 2525.94 s | dt: 27.742 s
Step: Lithiate | Increment: 219 | Iterations: 4
Simulation Time: 2553.68 s | dt: 27.742 s
Step: Lithiate | Increment: 220 | Iterations: 4
Simulation Time: 2581.42 s | dt: 27.742 s
Step: Lithiate | Increment: 221 | Iterations: 4
Simulation Time: 2609.16 s | dt: 27.742 s
Step: Lithiate | Increment: 222 | Iterations: 4
Simulation Time: 2636.9 s | dt: 27.742 s
Step: Lithiate | Increment: 223 | Iterations: 4
Simulation Time: 2664.64 s | dt: 27.742 s
Step: Lithiate | Increment: 224 | Iterations: 4
Simulation Time: 2692.39 s | dt: 27.742 s
Step: Lithiate | Increment: 225 | Iterations: 4
Simulation Time: 2720.13 s | dt: 27.742 s
Step: Lithiate | Increment: 226 | Iterations: 4
Simulation Time: 2747.87 s | dt: 27.742 s
Step: Lithiate | Increment: 227 | Iterations: 5
Simulation Time: 2775.61 s | dt: 27.742 s
Step: Lithiate | Increment: 228 | Iterations: 5
Simulation Time: 2803.35 s | dt: 27.742 s
Step: Lithiate | Increment: 229 | Iterations: 5
Simulation Time: 2831.1 s | dt: 27.742 s
Step: Lithiate | Increment: 230 | Iterations: 6
Simulation Time: 2858.84 s | dt: 27.742 s
2025-02-13 16:19:16.699 ( 74.390s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
2025-02-13 16:19:17.624 ( 75.315s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 231 | Iterations: 6
Simulation Time: 2865.77 s | dt: 6.935 s
2025-02-13 16:19:18.846 ( 76.537s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
2025-02-13 16:19:20.102 ( 77.793s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 232 | Iterations: 6
Simulation Time: 2867.51 s | dt: 1.734 s
2025-02-13 16:19:21.322 ( 79.013s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
2025-02-13 16:19:22.721 ( 80.412s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 233 | Iterations: 5
Simulation Time: 2867.94 s | dt: 0.433 s
Step: Lithiate | Increment: 234 | Iterations: 6
Simulation Time: 2868.37 s | dt: 0.433 s
2025-02-13 16:19:24.272 ( 81.963s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 235 | Iterations: 5
Simulation Time: 2868.59 s | dt: 0.217 s
2025-02-13 16:19:25.448 ( 83.139s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 236 | Iterations: 5
Simulation Time: 2868.7 s | dt: 0.108 s
Step: Lithiate | Increment: 237 | Iterations: 5
Simulation Time: 2868.81 s | dt: 0.108 s
Step: Lithiate | Increment: 238 | Iterations: 5
Simulation Time: 2868.91 s | dt: 0.108 s
Step: Lithiate | Increment: 239 | Iterations: 6
Simulation Time: 2869.02 s | dt: 0.108 s
Step: Lithiate | Increment: 240 | Iterations: 6
Simulation Time: 2869.13 s | dt: 0.108 s
Step: Lithiate | Increment: 241 | Iterations: 6
Simulation Time: 2869.24 s | dt: 0.108 s
Step: Lithiate | Increment: 242 | Iterations: 5
Simulation Time: 2869.35 s | dt: 0.108 s
Step: Lithiate | Increment: 243 | Iterations: 5
Simulation Time: 2869.46 s | dt: 0.108 s
Step: Lithiate | Increment: 244 | Iterations: 5
Simulation Time: 2869.57 s | dt: 0.108 s
Step: Lithiate | Increment: 245 | Iterations: 5
Simulation Time: 2869.67 s | dt: 0.108 s
Step: Lithiate | Increment: 246 | Iterations: 5
Simulation Time: 2869.78 s | dt: 0.108 s
Step: Lithiate | Increment: 247 | Iterations: 5
Simulation Time: 2869.89 s | dt: 0.108 s
Step: Lithiate | Increment: 248 | Iterations: 5
Simulation Time: 2870.0 s | dt: 0.108 s
Step: Lithiate | Increment: 249 | Iterations: 5
Simulation Time: 2870.11 s | dt: 0.108 s
Step: Lithiate | Increment: 250 | Iterations: 5
Simulation Time: 2870.22 s | dt: 0.108 s
Step: Lithiate | Increment: 251 | Iterations: 5
Simulation Time: 2870.32 s | dt: 0.108 s
Step: Lithiate | Increment: 252 | Iterations: 5
Simulation Time: 2870.43 s | dt: 0.108 s
Step: Lithiate | Increment: 253 | Iterations: 5
Simulation Time: 2870.54 s | dt: 0.108 s
Step: Lithiate | Increment: 254 | Iterations: 5
Simulation Time: 2870.65 s | dt: 0.108 s
Step: Lithiate | Increment: 255 | Iterations: 5
Simulation Time: 2870.76 s | dt: 0.108 s
Step: Lithiate | Increment: 256 | Iterations: 5
Simulation Time: 2870.87 s | dt: 0.108 s
Step: Lithiate | Increment: 257 | Iterations: 5
Simulation Time: 2870.97 s | dt: 0.108 s
Step: Lithiate | Increment: 258 | Iterations: 5
Simulation Time: 2871.08 s | dt: 0.108 s
Step: Lithiate | Increment: 259 | Iterations: 5
Simulation Time: 2871.19 s | dt: 0.108 s
Step: Lithiate | Increment: 260 | Iterations: 5
Simulation Time: 2871.3 s | dt: 0.108 s
Step: Lithiate | Increment: 261 | Iterations: 5
Simulation Time: 2871.41 s | dt: 0.108 s
Step: Lithiate | Increment: 262 | Iterations: 5
Simulation Time: 2871.52 s | dt: 0.108 s
Step: Lithiate | Increment: 263 | Iterations: 5
Simulation Time: 2871.62 s | dt: 0.108 s
Step: Lithiate | Increment: 264 | Iterations: 5
Simulation Time: 2871.73 s | dt: 0.108 s
Step: Lithiate | Increment: 265 | Iterations: 5
Simulation Time: 2871.84 s | dt: 0.108 s
Step: Lithiate | Increment: 266 | Iterations: 5
Simulation Time: 2871.95 s | dt: 0.108 s
Step: Lithiate | Increment: 267 | Iterations: 5
Simulation Time: 2872.06 s | dt: 0.108 s
Step: Lithiate | Increment: 268 | Iterations: 5
Simulation Time: 2872.17 s | dt: 0.108 s
Step: Lithiate | Increment: 269 | Iterations: 5
Simulation Time: 2872.27 s | dt: 0.108 s
Step: Lithiate | Increment: 270 | Iterations: 5
Simulation Time: 2872.38 s | dt: 0.108 s
Step: Lithiate | Increment: 271 | Iterations: 5
Simulation Time: 2872.49 s | dt: 0.108 s
Step: Lithiate | Increment: 272 | Iterations: 5
Simulation Time: 2872.6 s | dt: 0.108 s
Step: Lithiate | Increment: 273 | Iterations: 5
Simulation Time: 2872.71 s | dt: 0.108 s
Step: Lithiate | Increment: 274 | Iterations: 5
Simulation Time: 2872.82 s | dt: 0.108 s
Step: Lithiate | Increment: 275 | Iterations: 5
Simulation Time: 2872.92 s | dt: 0.108 s
Step: Lithiate | Increment: 276 | Iterations: 5
Simulation Time: 2873.03 s | dt: 0.108 s
Step: Lithiate | Increment: 277 | Iterations: 5
Simulation Time: 2873.14 s | dt: 0.108 s
Step: Lithiate | Increment: 278 | Iterations: 5
Simulation Time: 2873.25 s | dt: 0.108 s
Step: Lithiate | Increment: 279 | Iterations: 5
Simulation Time: 2873.36 s | dt: 0.108 s
Step: Lithiate | Increment: 280 | Iterations: 5
Simulation Time: 2873.47 s | dt: 0.108 s
Step: Lithiate | Increment: 281 | Iterations: 5
Simulation Time: 2873.57 s | dt: 0.108 s
Step: Lithiate | Increment: 282 | Iterations: 5
Simulation Time: 2873.68 s | dt: 0.108 s
Step: Lithiate | Increment: 283 | Iterations: 5
Simulation Time: 2873.79 s | dt: 0.108 s
Step: Lithiate | Increment: 284 | Iterations: 5
Simulation Time: 2873.9 s | dt: 0.108 s
Step: Lithiate | Increment: 285 | Iterations: 5
Simulation Time: 2874.01 s | dt: 0.108 s
Step: Lithiate | Increment: 286 | Iterations: 4
Simulation Time: 2874.12 s | dt: 0.108 s
Step: Lithiate | Increment: 287 | Iterations: 5
Simulation Time: 2874.28 s | dt: 0.163 s
Step: Lithiate | Increment: 288 | Iterations: 5
Simulation Time: 2874.44 s | dt: 0.163 s
Step: Lithiate | Increment: 289 | Iterations: 5
Simulation Time: 2874.6 s | dt: 0.163 s
Step: Lithiate | Increment: 290 | Iterations: 5
Simulation Time: 2874.77 s | dt: 0.163 s
Step: Lithiate | Increment: 291 | Iterations: 5
Simulation Time: 2874.93 s | dt: 0.163 s
Step: Lithiate | Increment: 292 | Iterations: 5
Simulation Time: 2875.09 s | dt: 0.163 s
Step: Lithiate | Increment: 293 | Iterations: 5
Simulation Time: 2875.25 s | dt: 0.163 s
Step: Lithiate | Increment: 294 | Iterations: 5
Simulation Time: 2875.42 s | dt: 0.163 s
Step: Lithiate | Increment: 295 | Iterations: 5
Simulation Time: 2875.58 s | dt: 0.163 s
Step: Lithiate | Increment: 296 | Iterations: 5
Simulation Time: 2875.74 s | dt: 0.163 s
Step: Lithiate | Increment: 297 | Iterations: 5
Simulation Time: 2875.9 s | dt: 0.163 s
Step: Lithiate | Increment: 298 | Iterations: 5
Simulation Time: 2876.07 s | dt: 0.163 s
Step: Lithiate | Increment: 299 | Iterations: 5
Simulation Time: 2876.23 s | dt: 0.163 s
Step: Lithiate | Increment: 300 | Iterations: 5
Simulation Time: 2876.39 s | dt: 0.163 s
Step: Lithiate | Increment: 301 | Iterations: 4
Simulation Time: 2876.55 s | dt: 0.163 s
Step: Lithiate | Increment: 302 | Iterations: 5
Simulation Time: 2876.8 s | dt: 0.244 s
Step: Lithiate | Increment: 303 | Iterations: 5
Simulation Time: 2877.04 s | dt: 0.244 s
Step: Lithiate | Increment: 304 | Iterations: 5
Simulation Time: 2877.29 s | dt: 0.244 s
Step: Lithiate | Increment: 305 | Iterations: 5
Simulation Time: 2877.53 s | dt: 0.244 s
Step: Lithiate | Increment: 306 | Iterations: 5
Simulation Time: 2877.77 s | dt: 0.244 s
Step: Lithiate | Increment: 307 | Iterations: 5
Simulation Time: 2878.02 s | dt: 0.244 s
Step: Lithiate | Increment: 308 | Iterations: 5
Simulation Time: 2878.26 s | dt: 0.244 s
Step: Lithiate | Increment: 309 | Iterations: 5
Simulation Time: 2878.51 s | dt: 0.244 s
Step: Lithiate | Increment: 310 | Iterations: 5
Simulation Time: 2878.75 s | dt: 0.244 s
Step: Lithiate | Increment: 311 | Iterations: 5
Simulation Time: 2878.99 s | dt: 0.244 s
Step: Lithiate | Increment: 312 | Iterations: 5
Simulation Time: 2879.24 s | dt: 0.244 s
Step: Lithiate | Increment: 313 | Iterations: 5
Simulation Time: 2879.48 s | dt: 0.244 s
Step: Lithiate | Increment: 314 | Iterations: 5
Simulation Time: 2879.72 s | dt: 0.244 s
Step: Lithiate | Increment: 315 | Iterations: 5
Simulation Time: 2879.97 s | dt: 0.244 s
Step: Lithiate | Increment: 316 | Iterations: 5
Simulation Time: 2880.21 s | dt: 0.244 s
Step: Lithiate | Increment: 317 | Iterations: 5
Simulation Time: 2880.46 s | dt: 0.244 s
Step: Lithiate | Increment: 318 | Iterations: 5
Simulation Time: 2880.7 s | dt: 0.244 s
Step: Lithiate | Increment: 319 | Iterations: 5
Simulation Time: 2880.94 s | dt: 0.244 s
Step: Lithiate | Increment: 320 | Iterations: 5
Simulation Time: 2881.19 s | dt: 0.244 s
Step: Lithiate | Increment: 321 | Iterations: 5
Simulation Time: 2881.43 s | dt: 0.244 s
Step: Lithiate | Increment: 322 | Iterations: 5
Simulation Time: 2881.68 s | dt: 0.244 s
Step: Lithiate | Increment: 323 | Iterations: 5
Simulation Time: 2881.92 s | dt: 0.244 s
Step: Lithiate | Increment: 324 | Iterations: 5
Simulation Time: 2882.16 s | dt: 0.244 s
Step: Lithiate | Increment: 325 | Iterations: 6
Simulation Time: 2882.41 s | dt: 0.244 s
Step: Lithiate | Increment: 326 | Iterations: 6
Simulation Time: 2882.65 s | dt: 0.244 s
2025-02-13 16:19:50.268 ( 107.959s) [main ] NewtonSolver.cpp:270 WARN| Newton solver did not converge.
Step: Lithiate | Increment: 327 | Iterations: 6
Simulation Time: 2882.77 s | dt: 0.122 s
Step: Lithiate | Increment: 328 | Iterations: 5
Simulation Time: 2882.89 s | dt: 0.122 s
Step: Lithiate | Increment: 329 | Iterations: 5
Simulation Time: 2883.02 s | dt: 0.122 s
Step: Lithiate | Increment: 330 | Iterations: 5
Simulation Time: 2883.14 s | dt: 0.122 s
Step: Lithiate | Increment: 331 | Iterations: 4
Simulation Time: 2883.26 s | dt: 0.122 s
Step: Lithiate | Increment: 332 | Iterations: 5
Simulation Time: 2883.44 s | dt: 0.183 s
Step: Lithiate | Increment: 333 | Iterations: 4
Simulation Time: 2883.63 s | dt: 0.183 s
Step: Lithiate | Increment: 334 | Iterations: 4
Simulation Time: 2883.9 s | dt: 0.274 s
Step: Lithiate | Increment: 335 | Iterations: 4
Simulation Time: 2884.31 s | dt: 0.411 s
Step: Lithiate | Increment: 336 | Iterations: 4
Simulation Time: 2884.93 s | dt: 0.617 s
Step: Lithiate | Increment: 337 | Iterations: 4
Simulation Time: 2885.85 s | dt: 0.926 s
Step: Lithiate | Increment: 338 | Iterations: 4
Simulation Time: 2887.24 s | dt: 1.389 s
Step: Lithiate | Increment: 339 | Iterations: 4
Simulation Time: 2889.33 s | dt: 2.083 s
Step: Lithiate | Increment: 340 | Iterations: 4
Simulation Time: 2892.45 s | dt: 3.124 s
Step: Lithiate | Increment: 341 | Iterations: 4
Simulation Time: 2897.14 s | dt: 4.687 s
Step: Lithiate | Increment: 342 | Iterations: 4
Simulation Time: 2904.17 s | dt: 7.03 s
Step: Lithiate | Increment: 343 | Iterations: 4
Simulation Time: 2914.71 s | dt: 10.545 s
Step: Lithiate | Increment: 344 | Iterations: 4
Simulation Time: 2930.53 s | dt: 15.818 s
Step: Lithiate | Increment: 345 | Iterations: 4
Simulation Time: 2954.26 s | dt: 23.726 s
Step: Lithiate | Increment: 346 | Iterations: 4
Simulation Time: 2989.85 s | dt: 35.59 s
Step: Lithiate | Increment: 347 | Iterations: 4
Simulation Time: 3025.44 s | dt: 35.59 s
Step: Lithiate | Increment: 348 | Iterations: 4
Simulation Time: 3061.03 s | dt: 35.59 s
Step: Lithiate | Increment: 349 | Iterations: 4
Simulation Time: 3096.62 s | dt: 35.59 s
Step: Lithiate | Increment: 350 | Iterations: 4
Simulation Time: 3132.2 s | dt: 35.59 s
Step: Lithiate | Increment: 351 | Iterations: 4
Simulation Time: 3167.79 s | dt: 35.59 s
Step: Lithiate | Increment: 352 | Iterations: 4
Simulation Time: 3203.38 s | dt: 35.59 s
Step: Lithiate | Increment: 353 | Iterations: 4
Simulation Time: 3238.97 s | dt: 35.59 s
Step: Lithiate | Increment: 354 | Iterations: 4
Simulation Time: 3274.56 s | dt: 35.59 s
Step: Lithiate | Increment: 355 | Iterations: 4
Simulation Time: 3310.15 s | dt: 35.59 s
Step: Lithiate | Increment: 356 | Iterations: 4
Simulation Time: 3345.74 s | dt: 35.59 s
Step: Lithiate | Increment: 357 | Iterations: 4
Simulation Time: 3381.33 s | dt: 35.59 s
Step: Lithiate | Increment: 358 | Iterations: 4
Simulation Time: 3416.92 s | dt: 35.59 s
Step: Lithiate | Increment: 359 | Iterations: 4
Simulation Time: 3452.51 s | dt: 35.59 s
Step: Lithiate | Increment: 360 | Iterations: 4
Simulation Time: 3488.1 s | dt: 35.59 s
Step: Lithiate | Increment: 361 | Iterations: 4
Simulation Time: 3523.69 s | dt: 35.59 s
Step: Lithiate | Increment: 362 | Iterations: 4
Simulation Time: 3559.28 s | dt: 35.59 s
Step: Lithiate | Increment: 363 | Iterations: 4
Simulation Time: 3594.87 s | dt: 35.59 s
Step: Lithiate | Increment: 364 | Iterations: 4
Simulation Time: 3630.46 s | dt: 35.59 s
-----------------------------------------
End computation
------------------------------------------
Elapsed real time: 0:01:49.875959
------------------------------------------
# set plot font to size 14
font = {'size' : 14}
plt.rc('font', **font)
# Only plot as far as we have time history data
ind = np.argmax(timeHist0)
# Get array of default plot colors
prop_cycle = plt.rcParams['axes.prop_cycle']
colors = prop_cycle.by_key()['color']
plt.figure()
plt.plot(timeHist2[0:ind], timeHist1[0:ind], linewidth=2.0,\
# color='r', marker='o', markersize=3)
color='b', label=r"$E/E_0 = 1.0$ ")
plt.axis('tight')
plt.xlabel(r"State of Charge")
plt.ylabel(r"$\hat\mu_{\text{ext}}= \bar\mu_{\text{ext}}/(R\theta)$")
plt.grid(linestyle="--", linewidth=0.5, color='b')
plt.xlim(-0.01,1.01)
plt.ylim(-2.5,2.5)
plt.legend()
fig = plt.gcf()
fig.set_size_inches(7,5)
plt.tight_layout()
plt.savefig("results/ellipse_lithiation_delithiation.png", dpi=600)

# plt.figure()
# plt.plot(timeHist1[0:ind], timeHist2[0:ind], linewidth=2.0,\
# # color='r', marker='o', markersize=3)
# color='b', label=r"$E/E_0 = 1.0$ ")
# plt.axis('tight')
# lt.xlabel(r"State of Charge")
# plt.ylabel(r"$\hat\mu= \mu/(R\theta)$")
# plt.grid(linestyle="--", linewidth=0.5, color='b')
# plt.ylim(0,1.01)
# plt.xlim(1060,1120)
# plt.legend()
# fig = plt.gcf()
# fig.set_size_inches(7,5)
# plt.tight_layout()
# plt.savefig("results/ellipse_lithiation_slow_ramp_detail_emod=1.0.png", dpi=600)