Simple stretch
Basic units: Length: mm Mass: kg Time: s Derived units: Force: milliNewtons Stress: kPa
Eric Stewart and Lallit Anand
ericstew@mit.edu and anand@mit.edu
Converted to FEniCSx by Jorge Nin jorgenin@mit.edu September 2023
import numpy as np
import dolfinx
from mpi4py import MPI
from petsc4py import PETSc
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 XDMFFile
import ufl
from ufl import (TestFunctions, TrialFunction, Identity, grad, det, div, dev, inv, tr, sqrt, conditional , gt, dx, inner, derivative, dot, ln, split)
from datetime import datetime
from dolfinx.plot import vtk_mesh
import pyvista
pyvista.set_jupyter_backend('client')
## Define temporal parameters
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)
length = 50.0 # mm
domain = mesh.create_box(MPI.COMM_WORLD,[[0.0,0.0,0.0],[length,length,length]],[10,10,6],mesh.CellType.tetrahedron)
plotter = pyvista.Plotter()
vtkdata = vtk_mesh(domain, domain.topology.dim)
grid = pyvista.UnstructuredGrid(*vtkdata)
actor = plotter.add_mesh(grid, show_edges=True)
plotter.show()
plotter.close()
Widget(value="<iframe src='http://localhost:51043/index.html?ui=P_0x2a4f0dcd0_0&reconnect=auto' style='width: …
def xBot(x):
return np.isclose(x[0], 0)
def xTop(x):
return np.isclose(x[0], length)
def yBot(x):
return np.isclose(x[1], 0)
def yTop(x):
return np.isclose(x[1], length)
def zBot(x):
return np.isclose(x[2], 0)
def loadFace(x):
return np.logical_and(np.logical_and(np.isclose(x[2],length) , np.less_equal(x[0],length/2)) , np.less_equal(x[1], length/2))
def zTop(x):
return np.logical_and(np.isclose(x[2], length),np.logical_not(loadFace(x)))
boundaries = [(2,xBot),(3,yBot),(4,zBot),(5,zTop),(6,loadFace)]
facet_indices, facet_markers = [], []
fdim = domain.topology.dim - 1
for (marker, locator) in boundaries:
facets = mesh.locate_entities(domain, fdim, locator)
facet_indices.append(facets)
facet_markers.append(np.full_like(facets, marker))
facet_indices = np.hstack(facet_indices).astype(np.int32)
facet_markers = np.hstack(facet_markers).astype(np.int32)
sorted_facets = np.argsort(facet_indices)
facet_tag = mesh.meshtags(domain, fdim, facet_indices[sorted_facets], facet_markers[sorted_facets])
ds = ufl.Measure('ds', domain=domain, subdomain_data=facet_tag,metadata={'quadrature_degree': 4})
n = ufl.FacetNormal(domain)
Arruda-Boyce Model
Gshear_0 = Constant(domain,PETSc.ScalarType(280.0)) # Ground state shear modulus
lambdaL = Constant(domain,PETSc.ScalarType(5.12)) # Locking stretch
Kbulk = Constant(domain,PETSc.ScalarType(1000.0*Gshear_0))
t = 0.0 # start time
Ttot = 30 # total simulation time
trac_tot = 1.5e3 # final pressure (kPa)
numSteps = 100
dt = Ttot/numSteps # fixed step size
def distRamp(t):
return trac_tot*t/Ttot
U2 = ufl.VectorElement("Lagrange", domain.ufl_cell(), 2) # For displacement
P1 = ufl.FiniteElement("Lagrange", domain.ufl_cell(), 1) # For pressure
TH = ufl.MixedElement([U2, P1]) # Taylor-Hood style mixed element
ME = FunctionSpace(domain, TH) # Total space for all DOFs
w = Function(ME)
u, p = split(w)
w_old = Function(ME)
u_old, p_old = split(w_old)
u_test, p_test = TestFunctions(ME)
dw = TrialFunction(ME)
def F_calc(u):
Id = Identity(3)
F = Id + grad(u)
return F
def lambdaBar_calc(u):
F = F_calc(u)
C = F.T*F
J = det(F)
Cdis = J**(-2/3)*C
I1 = tr(Cdis)
lambdaBar = sqrt(I1/3.0)
return lambdaBar
def zeta_calc(u):
lambdaBar = lambdaBar_calc(u)
# Use Pade approximation of Langevin inverse
z = lambdaBar/lambdaL
z = conditional(gt(z,0.95), 0.95, z) # Keep simulation from blowing up
beta = z*(3.0 - z**2.0)/(1.0 - z**2.0)
zeta = (lambdaL/(3*lambdaBar))*beta
return zeta
# Generalized shear modulus for Arruda-Boyce model
def Gshear_AB_calc(u):
zeta = zeta_calc(u)
Gshear = Gshear_0 * zeta
return Gshear
#---------------------------------------------
# Subroutine for calculating the Cauchy stress
#---------------------------------------------
def T_calc(u,p):
Id = Identity(3)
F = F_calc(u)
J = det(F)
B = F*F.T
Bdis = J**(-2/3)*B
Gshear = Gshear_AB_calc(u)
T = (1/J)* Gshear * dev(Bdis) - p * Id
return T
#----------------------------------------------
# Subroutine for calculating the Piola stress
#----------------------------------------------
def Piola_calc(u, p):
Id = Identity(3)
F = F_calc(u)
J = det(F)
#
T = T_calc(u,p)
#
Tmat = J * T * inv(F.T)
return Tmat
F = F_calc(u)
J = det(F)
lambdaBar = lambdaBar_calc(u)
# Piola stress
Tmat = Piola_calc(u, p)
dxs = dx(metadata={'quadrature_degree': 4})
# Residuals:
# Res_0: Balance of forces (test fxn: u)
# Res_1: Coupling pressure (test fxn: p)
Fcof = J*inv(F.T)
Time_cons = Constant(domain,PETSc.ScalarType(distRamp(0)))
traction = - Time_cons*dot(Fcof,n)
# The weak form for the equilibrium equation. No body force
Res_0 = inner(Tmat , grad(u_test) )*dxs- dot(traction, u_test)*ds(6)
# The weak form for the pressure
fac_p = ln(J)/J
#
Res_1 = dot( (p/Kbulk + fac_p), p_test)*dxs
# Total weak form
Res = Res_0 + Res_1
# Automatic differentiation tangent:
a = derivative(Res, w, dw)
#U0, submap = ME.sub(0).sub(1).collapse()
#fixed_displacement = fem.Function(U0)
#fixed_displacement.interpolate(lambda x : np.full(x.shape[1], distRamp(Time_cons)))
dofs1 = fem.locate_dofs_topological(ME.sub(0).sub(0), facet_tag.dim, facet_tag.find(2))
dofs2 = fem.locate_dofs_topological(ME.sub(0).sub(1), facet_tag.dim, facet_tag.find(3))
dofs3 = fem.locate_dofs_topological(ME.sub(0).sub(0), facet_tag.dim, facet_tag.find(6))
dofs4 = fem.locate_dofs_topological(ME.sub(0).sub(1), facet_tag.dim, facet_tag.find(6))
dofs5 = fem.locate_dofs_topological(ME.sub(0).sub(2), facet_tag.dim, facet_tag.find(4))
bcs_1 = dirichletbc(0.0, dofs1,ME.sub(0).sub(0)) # u1 fix - xBot
bcs_2 = dirichletbc(0.0, dofs2, ME.sub(0).sub(1)) # u2 fix - yBot
bcs_3 = dirichletbc(0.0,dofs3 ,ME.sub(0).sub(0)) # u3 fix - zBot
bcs_4 = dirichletbc(0.0, dofs4, ME.sub(0).sub(1)) # u2 fix - yBot
bcx_5 = dirichletbc(0.0,dofs5 ,ME.sub(0).sub(2)) # u3 fix - zBot
#
bcs = [bcs_1, bcs_2, bcs_3, bcs_4,bcx_5]
#Setting up visualziation
import pyvista
import matplotlib
import cmasher as cmr
import os
if not os.path.exists("results"):
# Create a new directory because it does not exist
os.makedirs("results")
pyvista.set_jupyter_backend('client')
plotter = pyvista.Plotter()
plotter.open_gif("results/3D_cube_footing.gif")
V = fem.FunctionSpace(domain,U2) ## Difference
u_n = fem.Function(V)
u_ex = Expression(w.sub(0),V.element.interpolation_points())
u_n.interpolate(u_ex)
topology, cells, geometry = plot.vtk_mesh(u_n.function_space)
function_grid = pyvista.UnstructuredGrid(topology, cells, geometry)
values = np.zeros((geometry.shape[0], 3))
#u0 = w.sub(0).collapse()
values[:, :len(u_n)] = u_n.x.array.reshape(geometry.shape[0], len(u_n))
function_grid["u"] = values
function_grid.set_active_vectors("u")
# Warp mesh by deformation
warped = function_grid.warp_by_vector("u", factor=0)
warped.set_active_vectors("u")
cmap = cmr.get_sub_cmap(matplotlib.colormaps.get_cmap("jet"), 0.1, 0.9,N=6)
# Add mesh to plotter and visualize
actor = plotter.add_mesh(warped, show_edges=True, lighting=False, clim=[0, 67.5],cmap=cmap)
# Compute magnitude of displacement to visualize in GIF
Vs = fem.FunctionSpace(domain, ("Lagrange", 2))
magnitude = fem.Function(Vs)
us = fem.Expression(u_n[1], Vs.element.interpolation_points())
magnitude.interpolate(us)
warped["mag"] = magnitude.x.array
#plotter.view_xy()
#plotter.camera.position=[5,25,200]
#plotter.camera.focal_point=[5,40,0]
plotter.write_frame()
## Functions for visualization
U1 = ufl.VectorElement("Lagrange", domain.ufl_cell(), 1)
V2 = fem.FunctionSpace(domain, U1)#Vector function space
V1 = fem.FunctionSpace(domain, P1)#Scalar function space
u_r = Function(V2)
u_r.name = "disp"
p_r = Function(V1)
p_r.name = "p"
J_vis = Function(V1)
J_vis.name = "J"
J_expr = Expression(J,V1.element.interpolation_points())
lambdaBar_Vis = Function(V1)
lambdaBar_Vis.name = "lambdaBar"
lambdaBar_expr = Expression(lambdaBar,V1.element.interpolation_points())
P11 = Function(V1)
P11.name = "P11"
P11_expr = Expression(Tmat[0,0],V1.element.interpolation_points())
P22 = Function(V1)
P22.name = "P22"
P22_expr = Expression(Tmat[1,1],V1.element.interpolation_points())
P33 = Function(V1)
P33.name = "P33"
P33_expr = Expression(Tmat[2,2],V1.element.interpolation_points())
T = Tmat*F.T/J
T0 = T - (1/3)*tr(T)*Identity(3)
Mises = sqrt((3/2)*inner(T0, T0))
Mises_Vis= Function(V1,name="Mises")
Mises_expr = Expression(Mises,V1.element.interpolation_points())
def InterpAndSave(t,file):
u_r.interpolate(w.sub(0))
p_r.interpolate(w.sub(1))
J_vis.interpolate(J_expr)
P11.interpolate(P11_expr)
P22.interpolate(P22_expr)
P33.interpolate(P33_expr)
lambdaBar_Vis.interpolate(lambdaBar_expr)
Mises_Vis.interpolate(Mises_expr)
file.write_function(u_r,t)
file.write_function(p_r,t)
file.write_function(J_vis,t)
file.write_function(P11,t)
file.write_function(P22,t)
file.write_function(P33,t)
file.write_function(lambdaBar_Vis,t)
file.write_function(Mises_Vis,t)
pointForStress = [0,0,length]
bb_tree = dolfinx.geometry.bb_tree(domain,domain.topology.dim)
cell_candidates = dolfinx.geometry.compute_collisions_points(bb_tree, pointForStress)
colliding_cells = dolfinx.geometry.compute_colliding_cells(domain, cell_candidates, pointForStress)
startTime = datetime.now()
#"-mcpu=apple-m1"
jit_options ={"cffi_extra_compile_args":["-O3","-ffast-math","-mcpu=apple-m1"]}
problem = NonlinearProblem(Res, w, bcs, a,jit_options=jit_options)
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
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"
#opts[f"{option_prefix}pc_hypre_type"] = "boomeramg"
opts[f"{option_prefix}pc_factor_mat_solver_type"] = "mumps"
opts[f"{option_prefix}ksp_max_it"] = 30
ksp.setFromOptions()
totSteps = numSteps+1
timeHist0 = np.zeros(shape=[totSteps])
step = "Strech"
t = 0
ii = 0
xdmf = XDMFFile(domain.comm, "results/3D_cube_footing.xdmf", "w")
xdmf.write_mesh(domain)
InterpAndSave(t,xdmf)
print("------------------------------------")
print("Simulation Start")
print("------------------------------------")
# Store start time
startTime = datetime.now()
while (round(t + dt, 9) <= Ttot):
# increment time
t += dt
# increment counter
ii += 1
# update time variables in time-dependent BCs
Time_cons.value = distRamp(t)
# Solve the problem
try:
(iter, converged) = solver.solve(w)
except: # Break the loop if solver fails
print( "Ended Early")
break
w.x.scatter_forward()
# Write output to *.xdmf file
#writeResults(t)
#print(u0.x.array-w.x.array[dofs])
#Visualizing GIF
u_n.interpolate(u_ex)
function_grid["u"][:, :len(u_n)] = u_n.x.array.reshape(geometry.shape[0], len(u_n))
magnitude.interpolate(us)
warped.set_active_scalars("mag")
warped_n = function_grid.warp_by_vector(factor=1)
plotter.update_coordinates(warped_n.points.copy(), render=False)
plotter.update_scalars(magnitude.x.array,render = False)
plotter.write_frame()
# Update DOFs for next step
w_old.x.array[:] = w.x.array
#SAVING RESULT
InterpAndSave(t,xdmf)
# Store displacement at a particular point at this time
timeHist0[ii] = w.sub(0).sub(2).eval([0, 0, length],colliding_cells[0])[0] # time history of displacement
#
# Print progress of calculation
if ii%1 == 0:
now = datetime.now()
current_time = now.strftime("%H:%M:%S")
print("Step: {} | Simulation Time: {} s, Wallclock Time: {}".\
format(step, round(t,4), current_time))
print("Iterations: {}".format(iter))
print()
plotter.close()
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("------------------------------------------")
xdmf.close()
------------------------------------ Simulation Start ------------------------------------ Step: Strech | Simulation Time: 0.3 s, Wallclock Time: 20:20:40 Iterations: 4 Step: Strech | Simulation Time: 0.6 s, Wallclock Time: 20:20:42 Iterations: 4 Step: Strech | Simulation Time: 0.9 s, Wallclock Time: 20:20:43 Iterations: 4 Step: Strech | Simulation Time: 1.2 s, Wallclock Time: 20:20:45 Iterations: 4 Step: Strech | Simulation Time: 1.5 s, Wallclock Time: 20:20:46 Iterations: 4 Step: Strech | Simulation Time: 1.8 s, Wallclock Time: 20:20:48 Iterations: 4 Step: Strech | Simulation Time: 2.1 s, Wallclock Time: 20:20:49 Iterations: 4 Step: Strech | Simulation Time: 2.4 s, Wallclock Time: 20:20:51 Iterations: 4 Step: Strech | Simulation Time: 2.7 s, Wallclock Time: 20:20:52 Iterations: 4 Step: Strech | Simulation Time: 3.0 s, Wallclock Time: 20:20:54 Iterations: 4 Step: Strech | Simulation Time: 3.3 s, Wallclock Time: 20:20:55 Iterations: 4 Step: Strech | Simulation Time: 3.6 s, Wallclock Time: 20:20:57 Iterations: 4 Step: Strech | Simulation Time: 3.9 s, Wallclock Time: 20:20:58 Iterations: 4 Step: Strech | Simulation Time: 4.2 s, Wallclock Time: 20:21:00 Iterations: 4 Step: Strech | Simulation Time: 4.5 s, Wallclock Time: 20:21:01 Iterations: 4 Step: Strech | Simulation Time: 4.8 s, Wallclock Time: 20:21:03 Iterations: 4 Step: Strech | Simulation Time: 5.1 s, Wallclock Time: 20:21:04 Iterations: 4 Step: Strech | Simulation Time: 5.4 s, Wallclock Time: 20:21:05 Iterations: 4 Step: Strech | Simulation Time: 5.7 s, Wallclock Time: 20:21:07 Iterations: 4 Step: Strech | Simulation Time: 6.0 s, Wallclock Time: 20:21:08 Iterations: 4 Step: Strech | Simulation Time: 6.3 s, Wallclock Time: 20:21:10 Iterations: 4 Step: Strech | Simulation Time: 6.6 s, Wallclock Time: 20:21:11 Iterations: 4 Step: Strech | Simulation Time: 6.9 s, Wallclock Time: 20:21:13 Iterations: 4 Step: Strech | Simulation Time: 7.2 s, Wallclock Time: 20:21:14 Iterations: 4 Step: Strech | Simulation Time: 7.5 s, Wallclock Time: 20:21:16 Iterations: 4 Step: Strech | Simulation Time: 7.8 s, Wallclock Time: 20:21:17 Iterations: 4 Step: Strech | Simulation Time: 8.1 s, Wallclock Time: 20:21:19 Iterations: 4 Step: Strech | Simulation Time: 8.4 s, Wallclock Time: 20:21:20 Iterations: 4 Step: Strech | Simulation Time: 8.7 s, Wallclock Time: 20:21:22 Iterations: 4 Step: Strech | Simulation Time: 9.0 s, Wallclock Time: 20:21:23 Iterations: 4 Step: Strech | Simulation Time: 9.3 s, Wallclock Time: 20:21:25 Iterations: 4 Step: Strech | Simulation Time: 9.6 s, Wallclock Time: 20:21:26 Iterations: 4 Step: Strech | Simulation Time: 9.9 s, Wallclock Time: 20:21:28 Iterations: 4 Step: Strech | Simulation Time: 10.2 s, Wallclock Time: 20:21:29 Iterations: 4 Step: Strech | Simulation Time: 10.5 s, Wallclock Time: 20:21:31 Iterations: 4 Step: Strech | Simulation Time: 10.8 s, Wallclock Time: 20:21:32 Iterations: 4 Step: Strech | Simulation Time: 11.1 s, Wallclock Time: 20:21:34 Iterations: 4 Step: Strech | Simulation Time: 11.4 s, Wallclock Time: 20:21:35 Iterations: 4 Step: Strech | Simulation Time: 11.7 s, Wallclock Time: 20:21:37 Iterations: 4 Step: Strech | Simulation Time: 12.0 s, Wallclock Time: 20:21:38 Iterations: 4 Step: Strech | Simulation Time: 12.3 s, Wallclock Time: 20:21:40 Iterations: 4 Step: Strech | Simulation Time: 12.6 s, Wallclock Time: 20:21:41 Iterations: 4 Step: Strech | Simulation Time: 12.9 s, Wallclock Time: 20:21:43 Iterations: 4 Step: Strech | Simulation Time: 13.2 s, Wallclock Time: 20:21:44 Iterations: 4 Step: Strech | Simulation Time: 13.5 s, Wallclock Time: 20:21:46 Iterations: 4 Step: Strech | Simulation Time: 13.8 s, Wallclock Time: 20:21:47 Iterations: 4 Step: Strech | Simulation Time: 14.1 s, Wallclock Time: 20:21:49 Iterations: 4 Step: Strech | Simulation Time: 14.4 s, Wallclock Time: 20:21:50 Iterations: 4 Step: Strech | Simulation Time: 14.7 s, Wallclock Time: 20:21:52 Iterations: 4 Step: Strech | Simulation Time: 15.0 s, Wallclock Time: 20:21:53 Iterations: 4 Step: Strech | Simulation Time: 15.3 s, Wallclock Time: 20:21:55 Iterations: 4 Step: Strech | Simulation Time: 15.6 s, Wallclock Time: 20:21:56 Iterations: 4 Step: Strech | Simulation Time: 15.9 s, Wallclock Time: 20:21:58 Iterations: 4 Step: Strech | Simulation Time: 16.2 s, Wallclock Time: 20:21:59 Iterations: 4 Step: Strech | Simulation Time: 16.5 s, Wallclock Time: 20:22:01 Iterations: 4 Step: Strech | Simulation Time: 16.8 s, Wallclock Time: 20:22:03 Iterations: 4 Step: Strech | Simulation Time: 17.1 s, Wallclock Time: 20:22:04 Iterations: 4 Step: Strech | Simulation Time: 17.4 s, Wallclock Time: 20:22:06 Iterations: 4 Step: Strech | Simulation Time: 17.7 s, Wallclock Time: 20:22:07 Iterations: 4 Step: Strech | Simulation Time: 18.0 s, Wallclock Time: 20:22:09 Iterations: 4 Step: Strech | Simulation Time: 18.3 s, Wallclock Time: 20:22:10 Iterations: 4 Step: Strech | Simulation Time: 18.6 s, Wallclock Time: 20:22:12 Iterations: 4 Step: Strech | Simulation Time: 18.9 s, Wallclock Time: 20:22:13 Iterations: 4 Step: Strech | Simulation Time: 19.2 s, Wallclock Time: 20:22:15 Iterations: 4 Step: Strech | Simulation Time: 19.5 s, Wallclock Time: 20:22:16 Iterations: 4 Step: Strech | Simulation Time: 19.8 s, Wallclock Time: 20:22:18 Iterations: 4 Step: Strech | Simulation Time: 20.1 s, Wallclock Time: 20:22:19 Iterations: 4 Step: Strech | Simulation Time: 20.4 s, Wallclock Time: 20:22:21 Iterations: 4 Step: Strech | Simulation Time: 20.7 s, Wallclock Time: 20:22:22 Iterations: 4 Step: Strech | Simulation Time: 21.0 s, Wallclock Time: 20:22:24 Iterations: 4 Step: Strech | Simulation Time: 21.3 s, Wallclock Time: 20:22:26 Iterations: 4 Step: Strech | Simulation Time: 21.6 s, Wallclock Time: 20:22:27 Iterations: 4 Step: Strech | Simulation Time: 21.9 s, Wallclock Time: 20:22:29 Iterations: 4 Step: Strech | Simulation Time: 22.2 s, Wallclock Time: 20:22:30 Iterations: 4 Step: Strech | Simulation Time: 22.5 s, Wallclock Time: 20:22:32 Iterations: 4 Step: Strech | Simulation Time: 22.8 s, Wallclock Time: 20:22:33 Iterations: 4 Step: Strech | Simulation Time: 23.1 s, Wallclock Time: 20:22:35 Iterations: 4 Step: Strech | Simulation Time: 23.4 s, Wallclock Time: 20:22:36 Iterations: 4 Step: Strech | Simulation Time: 23.7 s, Wallclock Time: 20:22:38 Iterations: 4 Step: Strech | Simulation Time: 24.0 s, Wallclock Time: 20:22:39 Iterations: 4 Step: Strech | Simulation Time: 24.3 s, Wallclock Time: 20:22:41 Iterations: 4 Step: Strech | Simulation Time: 24.6 s, Wallclock Time: 20:22:42 Iterations: 4 Step: Strech | Simulation Time: 24.9 s, Wallclock Time: 20:22:44 Iterations: 4 Step: Strech | Simulation Time: 25.2 s, Wallclock Time: 20:22:46 Iterations: 4 Step: Strech | Simulation Time: 25.5 s, Wallclock Time: 20:22:47 Iterations: 4 Step: Strech | Simulation Time: 25.8 s, Wallclock Time: 20:22:49 Iterations: 4 Step: Strech | Simulation Time: 26.1 s, Wallclock Time: 20:22:50 Iterations: 4 Step: Strech | Simulation Time: 26.4 s, Wallclock Time: 20:22:52 Iterations: 4 Step: Strech | Simulation Time: 26.7 s, Wallclock Time: 20:22:53 Iterations: 4 Step: Strech | Simulation Time: 27.0 s, Wallclock Time: 20:22:55 Iterations: 4 Step: Strech | Simulation Time: 27.3 s, Wallclock Time: 20:22:56 Iterations: 4 Step: Strech | Simulation Time: 27.6 s, Wallclock Time: 20:22:58 Iterations: 4 Step: Strech | Simulation Time: 27.9 s, Wallclock Time: 20:22:59 Iterations: 4 Step: Strech | Simulation Time: 28.2 s, Wallclock Time: 20:23:01 Iterations: 4 Step: Strech | Simulation Time: 28.5 s, Wallclock Time: 20:23:02 Iterations: 4 Step: Strech | Simulation Time: 28.8 s, Wallclock Time: 20:23:04 Iterations: 4 Step: Strech | Simulation Time: 29.1 s, Wallclock Time: 20:23:06 Iterations: 4 Step: Strech | Simulation Time: 29.4 s, Wallclock Time: 20:23:07 Iterations: 4 Step: Strech | Simulation Time: 29.7 s, Wallclock Time: 20:23:09 Iterations: 4 Step: Strech | Simulation Time: 30.0 s, Wallclock Time: 20:23:10 Iterations: 4 ----------------------------------------- End computation ------------------------------------------ Elapsed real time: 0:02:35.500574 ------------------------------------------
from IPython.display import Image
display(Image(data=open("results/3D_cube_footing.gif",'rb').read(), format='png'))