Thermal de-swelling of a square#
Swelling-deswelling-reswelling of a square domain
This is a two-dimensional plane-strain simulation
Degrees of freedom#
Displacement: u
pressure: p
chemical potential: mu
concentration: c
temperature: theta
Units#
Length: mm
Mass: tonne (1000kg)
Time: s
Mass density: tonne/mm^3
Force: N
Stress: MPa
Energy: mJ
Temperature: K
Amount of substance: mol
Species concentration: mol/mm^3
Chemical potential: mJ/mol
Molar volume: mm^3/mol
Species diffusivity: mm^2/s
Thermal expansion coefficient: #/K
Specific heat: mJ/(mm^3 K)
Thermal conductivity: mW/(mm K)
Boltzmann Constant: 1.38E-20 mJ/K
Gas constant: 8.3145E3 mJ/(mol 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 ,\
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#
# Create square mesh
L0 = 2.5 # mm, half-length of beam
# domain = mesh.create_rectangle(MPI.COMM_WORLD, [[0.0,0.0], [L0,L0]],\
# [15,15], mesh.CellType.quadrilateral)
domain = mesh.create_rectangle(MPI.COMM_WORLD, [[0.0,0.0], [L0,L0]],\
[15,15], mesh.CellType.triangle, diagonal=mesh.DiagonalType.crossed)
x = ufl.SpatialCoordinate(domain)
Identify the boundaries of the domain
# Identify the boundaries of the rectangle mesh
#
def xBot(x):
return np.isclose(x[0], 0)
def xTop(x):
return np.isclose(x[0], L0)
def yBot(x):
return np.isclose(x[1], 0)
def yTop(x):
return np.isclose(x[1], L0)
# Mark the sub-domains
boundaries = [(1, xBot),(2,xTop),(3,yBot),(4,yTop)]
# Build collections of facets on each subdomain and mark them appropriately.
facet_indices, facet_markers = [], [] # initalize empty collections of indices and markers.
fdim = domain.topology.dim - 1 # geometric dimension of the facet (mesh dimension - 1)
for (marker, locator) in boundaries:
facets = mesh.locate_entities(domain, fdim, locator) # an array of all the facets in a
# given subdomain ("locator")
facet_indices.append(facets) # add these facets to the collection.
facet_markers.append(np.full_like(facets, marker)) # mark them with the appropriate index.
# Format the facet indices and markers as required for use in dolfinx.
facet_indices = np.hstack(facet_indices).astype(np.int32)
facet_markers = np.hstack(facet_markers).astype(np.int32)
sorted_facets = np.argsort(facet_indices)
#
# Add these marked facets as "mesh tags" for later use in BCs.
facet_tags = mesh.meshtags(domain, fdim, facet_indices[sorted_facets], facet_markers[sorted_facets])
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/square_mesh.png")
from IPython.display import Image
Image(filename='results/square_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':2})
# Define the volume integration measure "dx"
# also specify the number of volume quadrature points.
dx = ufl.Measure('dx', domain=domain, metadata={'quadrature_degree': 2})
# 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#
# Set the locking stretch to a large number to model a Neo-Hookean material
#
k_B = Constant(domain,PETSc.ScalarType(1.38E-20)) # Boltzmann's constant
R_gas = Constant(domain,PETSc.ScalarType(8.3145E3)) # Gas constant
theta0 = Constant(domain,PETSc.ScalarType(298)) # Initial temperature
Gshear_0 = Constant(domain,PETSc.ScalarType(0.3)) # Ground sate shear modulus
N_R = Constant(domain,PETSc.ScalarType(Gshear_0/(k_B*theta0))) # Number polymer chains per unit ref. volume
lambdaL = Constant(domain,PETSc.ScalarType(5.2)) # Locking stretch
Kbulk = Constant(domain,PETSc.ScalarType(1000.0*Gshear_0)) # Bulk modulus
#
Omega = Constant(domain,PETSc.ScalarType(1.0E5)) # Molar volume of fluid
D = Constant(domain,PETSc.ScalarType(5.0E-3)) # Diffusivity
chi_L = Constant(domain,PETSc.ScalarType(0.1)) # Flory-Huggins mixing parameter low value
chi_H = Constant(domain,PETSc.ScalarType(0.7)) # Flory-Huggins mixing parameter high value
theta_T = Constant(domain,PETSc.ScalarType(307)) # Transition temperature
Delta_T = Constant(domain,PETSc.ScalarType(5.0)) # Transition temperature width
alpha = Constant(domain,PETSc.ScalarType(70.0E-6)) # Coefficient of thermal expansion
c_v = Constant(domain,PETSc.ScalarType(4.18)) # Specific heat
k_therm = Constant(domain,PETSc.ScalarType(0.53)) # Thermal conductivity
#
phi0 = Constant(domain,PETSc.ScalarType(0.999)) # Initial polymer volume fraction
mu0 = Constant(domain,PETSc.ScalarType(ln(1.0-phi0) + phi0 )) # Initial chemical potential
c0 = Constant(domain,PETSc.ScalarType((1/phi0) - 1)) # Initial concentration
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) # For pressure, chemical potential,
# species concentration, and temperature
#
TH = mixed_element([U2, P1, P1, P1, P1]) # Taylor-Hood style mixed element
ME = functionspace(domain, TH) # Total space for all DOFs
# Define actual functions with the required DOFs
w = Function(ME)
u, p, mu, c, theta = split(w) # displacement u, pressure p, chemical potential mu,
# concentration c, and temperature theta
# A copy of functions to store values in the previous step for time-stepping
w_old = Function(ME)
u_old, p_old, mu_old, c_old, theta_old = split(w_old)
# Define test functions
u_test, p_test, mu_test, c_test, theta_test = TestFunctions(ME)
# Define trial functions needed for automatic differentiation
dw = TrialFunction(ME)
Initial conditions#
The initial conditions for \(\mathbf{u}\) and \(p\) are zero everywhere.
These are imposed automatically, since we have not specified any non-zero initial conditions.
We do, however, need to impose the uniform initial conditions for \(\mu=\mu_0\) and \(\hat{c} = \hat{c}_0\) which correspond to \(\phi_0 = 0.999\), as well as the initial temperature \(\vartheta=\vartheta_0\). This is done below.
# Assign initial normalized chemical potential mu0 to the domain
w.sub(2).interpolate(lambda x: np.full((x.shape[1],), mu0))
w_old.sub(2).interpolate(lambda x: np.full((x.shape[1],), mu0))
# Assign initial value of normalized concentration c0 to the domain
w.sub(3).interpolate(lambda x: np.full((x.shape[1],), c0))
w_old.sub(3).interpolate(lambda x: np.full((x.shape[1],), c0))
# Assign initial value of tempearature theta0 to the domain
w.sub(4).interpolate(lambda x: np.full((x.shape[1],), theta0))
w_old.sub(4).interpolate(lambda x: np.full((x.shape[1],), theta0))
Subroutines for kinematics and constitutive equations#
# Special gradient operators for plane strain
#
# Gradient of vector field u
def pe_grad_vector(u):
grad_u = grad(u)
pe_grad_u = ufl.as_tensor([ [grad_u[0,0], grad_u[0,1], 0.0],
[grad_u[1,0], grad_u[1,1], 0.0],
[ 0.0, 0.0, 0.0] ])
return pe_grad_u
# Gradient of scalar field y
# (just need an extra zero for dimensions to work out)
def pe_grad_scalar(y):
grad_y = grad(y)
pe_grad_y = ufl.as_vector([grad_y[0], grad_y[1], 0.0])
return pe_grad_y
# Plane strain deformation gradient
def F_pe_calc(u):
dim = len(u) # dimension of problem (2)
Id = Identity(dim) # 2D Identity tensor
F = Id + grad(u) # 2D Deformation gradient
F_pe = ufl.as_tensor([ [F[0,0], F[0,1], 0.0],
[F[1,0], F[1,1], 0.0],
[ 0.0, 0.0, 1.0] ]) # Full plane strain F
return F_pe
# Elastic deformation gradient Fe
def Fe_calc(u,c):
F = F_pe_calc(u)
J = det(F)
#
Js = 1.0 + c
Fs = Js**(1/3)*Identity(3)
#
Fe = F*inv(Fs)
return Fe
def lambdaBar_calc(u):
F = F_pe_calc(u)
C = F.T*F
I1 = tr(C)
lambdaBar = sqrt(I1/3.0)
return lambdaBar
# Calculate zeta
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
# Calculate zeta0
def zeta0_calc():
# Use Pade approximation of Langevin inverse (A. Cohen, 1991)
z = 1/lambdaL
z = conditional(gt(z,0.95), 0.95, z) # Keep from blowing up
beta0 = z*(3.0 - z**2.0)/(1.0 - z**2.0)
zeta0 = (lambdaL/3)*beta0
return zeta0
# chi-parameter
def chi_calc(theta):
chi = 0.5*(chi_L + chi_H)- 0.5*(chi_L - chi_H)* tanh((theta-theta_T)/Delta_T)
return chi
# Stress-temperature modulus
def MH1_calc(u,c):
Id = Identity(3)
Fe = Fe_calc(u,c)
Je = det(Fe)
Ce = Fe.T*Fe
Js = 1 + c
zeta = zeta_calc(u)
zeta0 = zeta0_calc()
fac1 = N_R * k_B * ( zeta * Js**(2/3) * Id - zeta0 * inv(Ce) )
fac2 = (3*Kbulk*alpha/Je) * inv(Ce)
MH1 = 0.5*(fac1 +fac2)
return MH1
# Chemical potential-temperature modulus
def MH2_calc(u,c):
Id = Identity(3)
F = F_pe_calc(u)
C = F.T*F
Js = 1 + c
zeta = zeta_calc(u)
zeta0 = zeta0_calc()
chi = chi_calc(theta)
#
fac1 = R_gas*( ln(c/Js) + 1/Js + chi/Js**2 )
fac2 = (Omega*N_R*k_B/Js)*( zeta*tr(C)/3 - zeta0 )
MH2 = fac1+ fac2
return MH2
# Normalized Piola stress
def Piola_calc(u, p, theta):
F = F_pe_calc(u)
#
J = det(F)
#
zeta = zeta_calc(u)
#
zeta0 = zeta0_calc()
#
Gshear0 = N_R * k_B * theta
#
Piola = (zeta*F - zeta0*inv(F.T) ) - (J * p/Gshear0) * inv(F.T)
return Piola
# Species flux
def Flux_calc(u, mu, c, theta):
F = F_pe_calc(u)
#
Cinv = inv(F.T*F)
#
RT = R_gas * theta
#
Mob = (D*c)/(Omega*RT)*Cinv # Mobility tensor
#
Jmat = - RT* Mob * pe_grad_scalar(mu)
return Jmat
# Heat flux
def Heat_flux_calc(u, theta):
F = F_pe_calc(u)
J = det(F)
#
Cinv = inv(F.T*F)
#
Tcond = J * k_therm * Cinv # Thermal conductivity tensor
#
Qmat = - Tcond * pe_grad_scalar(theta)
return Qmat
Evaluate kinematics and constitutive relations#
# Kinematics
F = F_pe_calc(u)
J = det(F)
# lambdaBar
lambdaBar = lambdaBar_calc(u)
# Fe
Fe = Fe_calc(u,c)
Fe_old = Fe_calc(u_old,c_old)
# Je
Je = det(Fe)
# Ce
Ce = Fe.T * Fe
Ce_old = Fe_old.T * Fe_old
# Piola stress
Piola = Piola_calc(u, p, theta)
# Species flux
Jmat = Flux_calc(u, mu, c, theta)
# Heat flux
Qmat = Heat_flux_calc(u, theta)
# Heat-coupling terms
MH1 = MH1_calc(u,c)
MH2 = MH2_calc(u,c)
# RT
RT = R_gas*theta
# chi-parameter
chi = chi_calc(theta)
Weak forms#
# Residuals:
# Res_0: Balance of forces (test fxn: u)
# Res_1: Pressure (test fxn: p)
# Res_2: Balance of mass (test fxn: mu)
# Res_3: Concentration (test fxn:c)
# Res_4: Tempearture (test fxn:theta)
# The weak form for the equilibrium equation
Res_0 = inner(Piola , pe_grad_vector(u_test) )*dx
# The weak form for the pressure
fac_p = (ln(Je) - 3*alpha*(theta-theta0))
#
Res_1 = dot((p*Je/Kbulk + fac_p) , p_test)*dx
# The weak form for the mass balance of solvent
Res_2 = dot((c - c_old)/dk, mu_test)*dx \
- Omega*dot(Jmat , pe_grad_scalar(mu_test) )*dx
# The weak form for the concentration
fac = 1/(1+c)
#
fac1 = mu - ( ln(1.0-fac)+ fac + chi*fac*fac )
#
fac2 = -(Omega*Je/RT)*p
#
fac3 = - (Omega/RT)* ((1/2) * ((Je*p)**2.0)/Kbulk - 3*Kbulk*alpha*(theta-theta0)*ln(Je))
#
fac4 = fac1 + fac2 + fac3
#
Res_3 = dot(fac4, c_test)*dx
# The weak form for the the heat equation
tfac1 = theta * inner(MH1, (Ce-Ce_old))
tfac2 = theta * MH2*(c - c_old)/Omega
tfac3 = dk * RT * inner(Jmat, pe_grad_scalar(mu))
tfac4 = tfac1 + tfac2 - tfac3
#
Res_4 = dot( c_v*(theta - theta_old), theta_test)*dx \
- dk* dot(Qmat , pe_grad_scalar(theta_test) )*dx \
- dot(tfac4, theta_test)*dx
# Total weak form
Res = Res_0 + Res_1 + Res_2 + Res_3 + Res_4
# Automatic differentiation tangent:
a = derivative(Res, w, dw)
Set-up output files#
# Results file name
results_name = "thermogel_pe_square_swell"
# Function space for projection of results
U1 = element("DG", domain.basix_cell(), 1, shape=(2,)) # For displacement
P0 = element("DG", domain.basix_cell(), 1) # For pressure, chemical potential, and concentration
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"
p_vis = Function(V1)
p_vis.name = "p"
mu_vis = Function(V1)
mu_vis.name = "mu"
c_vis = Function(V1)
c_vis.name = "c"
theta_vis = Function(V1)
theta_vis.name = "theta"
# calculated fields to write to output file
phi = 1/(1+c)
phi_vis = Function(V1)
phi_vis.name = "phi"
phi_expr = Expression(phi,V1.element.interpolation_points())
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(Piola[0,0],V1.element.interpolation_points())
#
P22 = Function(V1)
P22.name = "P22"
P22_expr = Expression(Piola[1,1],V1.element.interpolation_points())
#
P33 = Function(V1)
P33.name = "P33"
P33_expr = Expression(Piola[2,2],V1.element.interpolation_points())
# Mises stress
T = Piola*F.T/J
T0 = T - (1/3)*tr(T)*Identity(3)
Cauchy_pressure = - (1/3)*tr(T)
Mises = sqrt((3/2)*inner(T0, T0))
#
Mises_vis = Function(V1,name="Mises")
Mises_expr = Expression(Mises,V1.element.interpolation_points())
#
Cauchy_pressure_vis = Function(V1,name="Cauchy_pressure")
Cauchy_pressure_expr = Expression(Cauchy_pressure,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, p_vis, mu_vis, c_vis, theta_vis, phi_vis, J_vis, P11, P22, P33,
lambdaBar_vis,Mises_vis,Cauchy_pressure_vis
],
engine="BP4",
)
def writeResults(t):
# Output field interpolation
u_vis.interpolate(w.sub(0))
p_vis.interpolate(w.sub(1))
mu_vis.interpolate(w.sub(2))
c_vis.interpolate(w.sub(3))
theta_vis.interpolate(w.sub(4))
phi_vis.interpolate(phi_expr)
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)
Cauchy_pressure_vis.interpolate(Cauchy_pressure_expr)
# Write output fields
file_results.write(t)
Infrastructure for pulling out time history data (displacement, force, etc.)#
# Code for visualizing concentration along a path
# (Based on Jeremy Bleyer's code for thermoleasticity)
#
# Collapse the concentration subspace:
V_C, _ = ME.sub(3).collapse()
#
# Identify concentration DOFs along desired plotting path which is yBot (index 3)
#
bottom_C_dofs = fem.locate_dofs_topological((ME.sub(3), V_C), facet_tags.dim, facet_tags.find(3))[1]
# Do not need this here
#
# # Identify point for reporting displacement
# pointForDisp = np.array([L0,L0/2,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
#
Analysis Step#
# Give the step a descriptive name
step = "Swell-deswell-reswell"
Boundary conditions#
# Boundary chemical potential history:
#
def muRamp(t):
return mu0*exp(-t/ttd)
# Constant for chemical potential history
#
mu_cons = Constant(domain,PETSc.ScalarType(muRamp(0)))
# t = 0.0 # initialization of time
# ttd = 300 # time for increasing mu
# tRamp = 3600*2 # time to increase or decrease tempearture
# time1 = 3600*6
# #
# time2 = 3600*8
# #
# time3 = 3600*14
# #
# time4 = 3600*16
# #
# time5 = 3600*20
# Boundary temperature history:
#
# Temperature increase deltaTheta =25K
#
deltaTheta = Constant(domain,PETSc.ScalarType(25))
#
def thetaRamp(t):
if t<= time1:
temperature = theta0
elif time1 < t <= time2:
temperature = theta0 + deltaTheta*(t-time1)/tRamp
elif time2 < t <= time3:
temperature = (theta0 + deltaTheta)
elif time3 < t <= time4 :
temperature = (theta0 + deltaTheta) - deltaTheta*(t-time3)/tRamp
else:
temperature = theta0
return temperature
# Constant for temperature history
theta_cons = Constant(domain,PETSc.ScalarType(thetaRamp(0)))
# Recall the sub-domains names and numbers
# boundaries = [(1, xBot),(2,xTop),(3,yBot),(4,yTop)]
# 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(1))
yBot_u2_dofs = fem.locate_dofs_topological(ME.sub(0).sub(1), facet_tags.dim, facet_tags.find(3))
#
xTop_mu_dofs = fem.locate_dofs_topological(ME.sub(2), facet_tags.dim, facet_tags.find(2))
yTop_mu_dofs = fem.locate_dofs_topological(ME.sub(2), facet_tags.dim, facet_tags.find(4))
#
xTop_theta_dofs = fem.locate_dofs_topological(ME.sub(4), facet_tags.dim, facet_tags.find(2))
yTop_theta_dofs = fem.locate_dofs_topological(ME.sub(4), facet_tags.dim, facet_tags.find(4))
# Build Dirichlet BCs
bcs_1 = dirichletbc(0.0, xBot_u1_dofs, ME.sub(0).sub(0)) # u1 fix - xBot
bcs_2 = dirichletbc(0.0, yBot_u2_dofs, ME.sub(0).sub(1)) # u2 fix - yBot
#
bcs_3 = dirichletbc(mu_cons, xTop_mu_dofs, ME.sub(2)) # mu_cons - xTop
bcs_4 = dirichletbc(mu_cons, yTop_mu_dofs, ME.sub(2)) # mu_cons - yTop
#
bcs_5 = dirichletbc(theta_cons, xTop_theta_dofs, ME.sub(4)) # theta_cons- xTop
bcs_6 = dirichletbc(theta_cons, yTop_theta_dofs, ME.sub(4)) # theta_cons - yTop
bcs = [bcs_1, bcs_2, bcs_3, bcs_4, bcs_5, bcs_6]
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
# 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])
# #
# timeHist3[0] = mu0 # Initialize the chemical potential
# Initialize a counter for reporting data
ii = 0
iii = 0
# Write initial state to file
writeResults(t=0.0)
# Initialize an array for storing results from a path:
#
xPath = V_C.tabulate_dof_coordinates()[bottom_C_dofs, 0] # x position of dofs
#
# Times for storing the profile data:
profileTimes = 3600*np.array([2, 6, 8, 14, 16, 20])
#
# Array for storing path data at specified times
#
C_res = np.zeros((len(xPath), len(profileTimes)))
Start calculation loop#
print("------------------------------------")
print("Simulation start")
print("------------------------------------")
# Store start time
startTime = datetime.now()
# Time-stepping solution procedure loop
while (round(t + dt, 9) <= Ttot):
# increment time
t += dt
# increment counter
ii += 1
# update time variables in time-dependent BCs
mu_cons.value = float(muRamp(t))
theta_cons.value = float(thetaRamp(t))
# Solve the problem
try:
(iter, converged) = solver.solve(w)
except: # Break the loop if solver fails
print("Ended Early")
break
# Collect results from MPI ghost processes
w.x.scatter_forward()
# Write output to file
writeResults(t)
# store profile data at desired times
if int(t) in profileTimes:
# write concentration profile data
C_out = w.sub(3).collapse()
C_out.name = "Concentration"
#
C_res[:, iii] = C_out.vector.array[bottom_C_dofs]
#
iii += 1
# Update DOFs for next step
w_old.x.array[:] = w.x.array
# Print progress of calculation
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 of {} s".\
format(round(t,4), Ttot))
print()
# 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: Swell-deswell-reswell | Increment: 1, Iterations: 6
Simulation Time: 100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 2, Iterations: 6
Simulation Time: 200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 3, Iterations: 6
Simulation Time: 300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 4, Iterations: 5
Simulation Time: 400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 5, Iterations: 5
Simulation Time: 500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 6, Iterations: 5
Simulation Time: 600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 7, Iterations: 5
Simulation Time: 700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 8, Iterations: 5
Simulation Time: 800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 9, Iterations: 5
Simulation Time: 900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 10, Iterations: 5
Simulation Time: 1000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 11, Iterations: 5
Simulation Time: 1100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 12, Iterations: 5
Simulation Time: 1200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 13, Iterations: 5
Simulation Time: 1300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 14, Iterations: 5
Simulation Time: 1400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 15, Iterations: 5
Simulation Time: 1500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 16, Iterations: 5
Simulation Time: 1600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 17, Iterations: 5
Simulation Time: 1700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 18, Iterations: 5
Simulation Time: 1800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 19, Iterations: 5
Simulation Time: 1900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 20, Iterations: 4
Simulation Time: 2000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 21, Iterations: 4
Simulation Time: 2100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 22, Iterations: 4
Simulation Time: 2200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 23, Iterations: 4
Simulation Time: 2300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 24, Iterations: 4
Simulation Time: 2400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 25, Iterations: 4
Simulation Time: 2500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 26, Iterations: 4
Simulation Time: 2600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 27, Iterations: 4
Simulation Time: 2700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 28, Iterations: 4
Simulation Time: 2800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 29, Iterations: 4
Simulation Time: 2900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 30, Iterations: 4
Simulation Time: 3000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 31, Iterations: 4
Simulation Time: 3100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 32, Iterations: 4
Simulation Time: 3200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 33, Iterations: 4
Simulation Time: 3300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 34, Iterations: 4
Simulation Time: 3400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 35, Iterations: 4
Simulation Time: 3500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 36, Iterations: 4
Simulation Time: 3600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 37, Iterations: 4
Simulation Time: 3700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 38, Iterations: 4
Simulation Time: 3800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 39, Iterations: 4
Simulation Time: 3900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 40, Iterations: 4
Simulation Time: 4000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 41, Iterations: 4
Simulation Time: 4100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 42, Iterations: 4
Simulation Time: 4200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 43, Iterations: 4
Simulation Time: 4300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 44, Iterations: 4
Simulation Time: 4400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 45, Iterations: 4
Simulation Time: 4500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 46, Iterations: 3
Simulation Time: 4600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 47, Iterations: 3
Simulation Time: 4700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 48, Iterations: 3
Simulation Time: 4800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 49, Iterations: 3
Simulation Time: 4900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 50, Iterations: 3
Simulation Time: 5000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 51, Iterations: 3
Simulation Time: 5100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 52, Iterations: 3
Simulation Time: 5200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 53, Iterations: 3
Simulation Time: 5300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 54, Iterations: 3
Simulation Time: 5400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 55, Iterations: 3
Simulation Time: 5500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 56, Iterations: 3
Simulation Time: 5600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 57, Iterations: 4
Simulation Time: 5700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 58, Iterations: 3
Simulation Time: 5800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 59, Iterations: 3
Simulation Time: 5900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 60, Iterations: 3
Simulation Time: 6000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 61, Iterations: 3
Simulation Time: 6100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 62, Iterations: 3
Simulation Time: 6200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 63, Iterations: 3
Simulation Time: 6300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 64, Iterations: 3
Simulation Time: 6400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 65, Iterations: 3
Simulation Time: 6500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 66, Iterations: 3
Simulation Time: 6600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 67, Iterations: 3
Simulation Time: 6700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 68, Iterations: 3
Simulation Time: 6800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 69, Iterations: 3
Simulation Time: 6900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 70, Iterations: 3
Simulation Time: 7000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 71, Iterations: 3
Simulation Time: 7100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 72, Iterations: 3
Simulation Time: 7200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 73, Iterations: 3
Simulation Time: 7300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 74, Iterations: 3
Simulation Time: 7400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 75, Iterations: 3
Simulation Time: 7500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 76, Iterations: 3
Simulation Time: 7600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 77, Iterations: 3
Simulation Time: 7700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 78, Iterations: 3
Simulation Time: 7800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 79, Iterations: 3
Simulation Time: 7900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 80, Iterations: 3
Simulation Time: 8000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 81, Iterations: 3
Simulation Time: 8100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 82, Iterations: 3
Simulation Time: 8200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 83, Iterations: 3
Simulation Time: 8300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 84, Iterations: 3
Simulation Time: 8400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 85, Iterations: 3
Simulation Time: 8500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 86, Iterations: 3
Simulation Time: 8600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 87, Iterations: 3
Simulation Time: 8700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 88, Iterations: 3
Simulation Time: 8800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 89, Iterations: 3
Simulation Time: 8900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 90, Iterations: 3
Simulation Time: 9000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 91, Iterations: 3
Simulation Time: 9100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 92, Iterations: 3
Simulation Time: 9200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 93, Iterations: 3
Simulation Time: 9300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 94, Iterations: 3
Simulation Time: 9400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 95, Iterations: 3
Simulation Time: 9500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 96, Iterations: 3
Simulation Time: 9600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 97, Iterations: 3
Simulation Time: 9700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 98, Iterations: 3
Simulation Time: 9800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 99, Iterations: 3
Simulation Time: 9900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 100, Iterations: 3
Simulation Time: 10000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 101, Iterations: 3
Simulation Time: 10100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 102, Iterations: 3
Simulation Time: 10200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 103, Iterations: 3
Simulation Time: 10300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 104, Iterations: 3
Simulation Time: 10400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 105, Iterations: 3
Simulation Time: 10500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 106, Iterations: 3
Simulation Time: 10600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 107, Iterations: 3
Simulation Time: 10700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 108, Iterations: 3
Simulation Time: 10800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 109, Iterations: 3
Simulation Time: 10900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 110, Iterations: 3
Simulation Time: 11000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 111, Iterations: 3
Simulation Time: 11100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 112, Iterations: 3
Simulation Time: 11200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 113, Iterations: 3
Simulation Time: 11300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 114, Iterations: 3
Simulation Time: 11400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 115, Iterations: 3
Simulation Time: 11500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 116, Iterations: 3
Simulation Time: 11600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 117, Iterations: 3
Simulation Time: 11700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 118, Iterations: 3
Simulation Time: 11800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 119, Iterations: 3
Simulation Time: 11900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 120, Iterations: 3
Simulation Time: 12000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 121, Iterations: 3
Simulation Time: 12100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 122, Iterations: 3
Simulation Time: 12200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 123, Iterations: 3
Simulation Time: 12300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 124, Iterations: 3
Simulation Time: 12400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 125, Iterations: 3
Simulation Time: 12500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 126, Iterations: 3
Simulation Time: 12600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 127, Iterations: 3
Simulation Time: 12700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 128, Iterations: 3
Simulation Time: 12800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 129, Iterations: 3
Simulation Time: 12900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 130, Iterations: 3
Simulation Time: 13000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 131, Iterations: 3
Simulation Time: 13100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 132, Iterations: 3
Simulation Time: 13200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 133, Iterations: 3
Simulation Time: 13300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 134, Iterations: 3
Simulation Time: 13400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 135, Iterations: 3
Simulation Time: 13500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 136, Iterations: 3
Simulation Time: 13600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 137, Iterations: 3
Simulation Time: 13700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 138, Iterations: 3
Simulation Time: 13800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 139, Iterations: 3
Simulation Time: 13900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 140, Iterations: 3
Simulation Time: 14000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 141, Iterations: 3
Simulation Time: 14100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 142, Iterations: 3
Simulation Time: 14200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 143, Iterations: 3
Simulation Time: 14300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 144, Iterations: 3
Simulation Time: 14400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 145, Iterations: 3
Simulation Time: 14500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 146, Iterations: 3
Simulation Time: 14600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 147, Iterations: 3
Simulation Time: 14700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 148, Iterations: 3
Simulation Time: 14800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 149, Iterations: 3
Simulation Time: 14900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 150, Iterations: 3
Simulation Time: 15000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 151, Iterations: 3
Simulation Time: 15100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 152, Iterations: 3
Simulation Time: 15200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 153, Iterations: 3
Simulation Time: 15300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 154, Iterations: 3
Simulation Time: 15400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 155, Iterations: 3
Simulation Time: 15500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 156, Iterations: 3
Simulation Time: 15600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 157, Iterations: 3
Simulation Time: 15700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 158, Iterations: 3
Simulation Time: 15800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 159, Iterations: 3
Simulation Time: 15900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 160, Iterations: 3
Simulation Time: 16000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 161, Iterations: 3
Simulation Time: 16100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 162, Iterations: 3
Simulation Time: 16200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 163, Iterations: 3
Simulation Time: 16300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 164, Iterations: 3
Simulation Time: 16400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 165, Iterations: 3
Simulation Time: 16500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 166, Iterations: 3
Simulation Time: 16600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 167, Iterations: 3
Simulation Time: 16700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 168, Iterations: 3
Simulation Time: 16800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 169, Iterations: 3
Simulation Time: 16900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 170, Iterations: 3
Simulation Time: 17000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 171, Iterations: 3
Simulation Time: 17100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 172, Iterations: 3
Simulation Time: 17200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 173, Iterations: 3
Simulation Time: 17300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 174, Iterations: 3
Simulation Time: 17400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 175, Iterations: 3
Simulation Time: 17500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 176, Iterations: 3
Simulation Time: 17600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 177, Iterations: 3
Simulation Time: 17700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 178, Iterations: 3
Simulation Time: 17800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 179, Iterations: 3
Simulation Time: 17900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 180, Iterations: 3
Simulation Time: 18000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 181, Iterations: 3
Simulation Time: 18100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 182, Iterations: 3
Simulation Time: 18200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 183, Iterations: 3
Simulation Time: 18300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 184, Iterations: 3
Simulation Time: 18400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 185, Iterations: 3
Simulation Time: 18500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 186, Iterations: 3
Simulation Time: 18600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 187, Iterations: 3
Simulation Time: 18700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 188, Iterations: 3
Simulation Time: 18800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 189, Iterations: 3
Simulation Time: 18900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 190, Iterations: 3
Simulation Time: 19000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 191, Iterations: 3
Simulation Time: 19100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 192, Iterations: 3
Simulation Time: 19200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 193, Iterations: 3
Simulation Time: 19300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 194, Iterations: 3
Simulation Time: 19400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 195, Iterations: 3
Simulation Time: 19500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 196, Iterations: 3
Simulation Time: 19600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 197, Iterations: 3
Simulation Time: 19700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 198, Iterations: 3
Simulation Time: 19800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 199, Iterations: 3
Simulation Time: 19900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 200, Iterations: 3
Simulation Time: 20000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 201, Iterations: 3
Simulation Time: 20100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 202, Iterations: 3
Simulation Time: 20200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 203, Iterations: 3
Simulation Time: 20300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 204, Iterations: 3
Simulation Time: 20400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 205, Iterations: 3
Simulation Time: 20500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 206, Iterations: 3
Simulation Time: 20600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 207, Iterations: 3
Simulation Time: 20700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 208, Iterations: 3
Simulation Time: 20800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 209, Iterations: 3
Simulation Time: 20900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 210, Iterations: 3
Simulation Time: 21000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 211, Iterations: 3
Simulation Time: 21100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 212, Iterations: 3
Simulation Time: 21200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 213, Iterations: 3
Simulation Time: 21300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 214, Iterations: 3
Simulation Time: 21400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 215, Iterations: 3
Simulation Time: 21500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 216, Iterations: 3
Simulation Time: 21600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 217, Iterations: 4
Simulation Time: 21700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 218, Iterations: 4
Simulation Time: 21800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 219, Iterations: 4
Simulation Time: 21900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 220, Iterations: 4
Simulation Time: 22000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 221, Iterations: 4
Simulation Time: 22100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 222, Iterations: 4
Simulation Time: 22200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 223, Iterations: 4
Simulation Time: 22300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 224, Iterations: 4
Simulation Time: 22400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 225, Iterations: 4
Simulation Time: 22500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 226, Iterations: 4
Simulation Time: 22600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 227, Iterations: 4
Simulation Time: 22700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 228, Iterations: 4
Simulation Time: 22800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 229, Iterations: 4
Simulation Time: 22900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 230, Iterations: 4
Simulation Time: 23000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 231, Iterations: 4
Simulation Time: 23100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 232, Iterations: 4
Simulation Time: 23200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 233, Iterations: 4
Simulation Time: 23300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 234, Iterations: 4
Simulation Time: 23400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 235, Iterations: 4
Simulation Time: 23500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 236, Iterations: 4
Simulation Time: 23600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 237, Iterations: 4
Simulation Time: 23700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 238, Iterations: 4
Simulation Time: 23800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 239, Iterations: 4
Simulation Time: 23900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 240, Iterations: 4
Simulation Time: 24000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 241, Iterations: 4
Simulation Time: 24100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 242, Iterations: 4
Simulation Time: 24200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 243, Iterations: 4
Simulation Time: 24300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 244, Iterations: 4
Simulation Time: 24400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 245, Iterations: 4
Simulation Time: 24500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 246, Iterations: 4
Simulation Time: 24600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 247, Iterations: 4
Simulation Time: 24700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 248, Iterations: 4
Simulation Time: 24800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 249, Iterations: 4
Simulation Time: 24900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 250, Iterations: 4
Simulation Time: 25000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 251, Iterations: 4
Simulation Time: 25100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 252, Iterations: 4
Simulation Time: 25200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 253, Iterations: 4
Simulation Time: 25300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 254, Iterations: 4
Simulation Time: 25400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 255, Iterations: 4
Simulation Time: 25500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 256, Iterations: 4
Simulation Time: 25600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 257, Iterations: 4
Simulation Time: 25700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 258, Iterations: 4
Simulation Time: 25800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 259, Iterations: 4
Simulation Time: 25900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 260, Iterations: 4
Simulation Time: 26000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 261, Iterations: 4
Simulation Time: 26100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 262, Iterations: 4
Simulation Time: 26200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 263, Iterations: 4
Simulation Time: 26300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 264, Iterations: 4
Simulation Time: 26400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 265, Iterations: 4
Simulation Time: 26500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 266, Iterations: 4
Simulation Time: 26600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 267, Iterations: 4
Simulation Time: 26700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 268, Iterations: 4
Simulation Time: 26800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 269, Iterations: 4
Simulation Time: 26900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 270, Iterations: 4
Simulation Time: 27000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 271, Iterations: 4
Simulation Time: 27100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 272, Iterations: 4
Simulation Time: 27200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 273, Iterations: 4
Simulation Time: 27300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 274, Iterations: 4
Simulation Time: 27400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 275, Iterations: 4
Simulation Time: 27500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 276, Iterations: 4
Simulation Time: 27600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 277, Iterations: 4
Simulation Time: 27700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 278, Iterations: 4
Simulation Time: 27800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 279, Iterations: 4
Simulation Time: 27900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 280, Iterations: 4
Simulation Time: 28000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 281, Iterations: 4
Simulation Time: 28100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 282, Iterations: 4
Simulation Time: 28200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 283, Iterations: 4
Simulation Time: 28300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 284, Iterations: 4
Simulation Time: 28400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 285, Iterations: 4
Simulation Time: 28500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 286, Iterations: 4
Simulation Time: 28600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 287, Iterations: 4
Simulation Time: 28700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 288, Iterations: 4
Simulation Time: 28800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 289, Iterations: 4
Simulation Time: 28900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 290, Iterations: 4
Simulation Time: 29000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 291, Iterations: 4
Simulation Time: 29100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 292, Iterations: 4
Simulation Time: 29200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 293, Iterations: 4
Simulation Time: 29300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 294, Iterations: 4
Simulation Time: 29400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 295, Iterations: 4
Simulation Time: 29500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 296, Iterations: 4
Simulation Time: 29600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 297, Iterations: 4
Simulation Time: 29700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 298, Iterations: 4
Simulation Time: 29800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 299, Iterations: 4
Simulation Time: 29900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 300, Iterations: 4
Simulation Time: 30000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 301, Iterations: 4
Simulation Time: 30100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 302, Iterations: 4
Simulation Time: 30200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 303, Iterations: 4
Simulation Time: 30300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 304, Iterations: 4
Simulation Time: 30400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 305, Iterations: 4
Simulation Time: 30500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 306, Iterations: 4
Simulation Time: 30600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 307, Iterations: 4
Simulation Time: 30700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 308, Iterations: 4
Simulation Time: 30800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 309, Iterations: 4
Simulation Time: 30900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 310, Iterations: 4
Simulation Time: 31000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 311, Iterations: 4
Simulation Time: 31100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 312, Iterations: 4
Simulation Time: 31200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 313, Iterations: 4
Simulation Time: 31300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 314, Iterations: 4
Simulation Time: 31400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 315, Iterations: 4
Simulation Time: 31500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 316, Iterations: 3
Simulation Time: 31600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 317, Iterations: 3
Simulation Time: 31700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 318, Iterations: 3
Simulation Time: 31800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 319, Iterations: 3
Simulation Time: 31900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 320, Iterations: 3
Simulation Time: 32000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 321, Iterations: 3
Simulation Time: 32100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 322, Iterations: 3
Simulation Time: 32200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 323, Iterations: 3
Simulation Time: 32300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 324, Iterations: 3
Simulation Time: 32400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 325, Iterations: 3
Simulation Time: 32500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 326, Iterations: 3
Simulation Time: 32600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 327, Iterations: 3
Simulation Time: 32700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 328, Iterations: 3
Simulation Time: 32800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 329, Iterations: 3
Simulation Time: 32900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 330, Iterations: 3
Simulation Time: 33000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 331, Iterations: 3
Simulation Time: 33100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 332, Iterations: 3
Simulation Time: 33200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 333, Iterations: 3
Simulation Time: 33300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 334, Iterations: 3
Simulation Time: 33400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 335, Iterations: 3
Simulation Time: 33500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 336, Iterations: 3
Simulation Time: 33600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 337, Iterations: 3
Simulation Time: 33700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 338, Iterations: 3
Simulation Time: 33800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 339, Iterations: 3
Simulation Time: 33900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 340, Iterations: 3
Simulation Time: 34000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 341, Iterations: 3
Simulation Time: 34100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 342, Iterations: 3
Simulation Time: 34200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 343, Iterations: 3
Simulation Time: 34300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 344, Iterations: 3
Simulation Time: 34400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 345, Iterations: 3
Simulation Time: 34500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 346, Iterations: 3
Simulation Time: 34600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 347, Iterations: 3
Simulation Time: 34700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 348, Iterations: 3
Simulation Time: 34800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 349, Iterations: 3
Simulation Time: 34900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 350, Iterations: 3
Simulation Time: 35000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 351, Iterations: 3
Simulation Time: 35100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 352, Iterations: 3
Simulation Time: 35200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 353, Iterations: 3
Simulation Time: 35300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 354, Iterations: 3
Simulation Time: 35400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 355, Iterations: 3
Simulation Time: 35500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 356, Iterations: 3
Simulation Time: 35600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 357, Iterations: 3
Simulation Time: 35700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 358, Iterations: 3
Simulation Time: 35800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 359, Iterations: 3
Simulation Time: 35900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 360, Iterations: 3
Simulation Time: 36000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 361, Iterations: 3
Simulation Time: 36100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 362, Iterations: 3
Simulation Time: 36200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 363, Iterations: 3
Simulation Time: 36300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 364, Iterations: 3
Simulation Time: 36400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 365, Iterations: 3
Simulation Time: 36500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 366, Iterations: 3
Simulation Time: 36600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 367, Iterations: 3
Simulation Time: 36700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 368, Iterations: 3
Simulation Time: 36800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 369, Iterations: 3
Simulation Time: 36900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 370, Iterations: 3
Simulation Time: 37000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 371, Iterations: 3
Simulation Time: 37100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 372, Iterations: 3
Simulation Time: 37200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 373, Iterations: 3
Simulation Time: 37300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 374, Iterations: 3
Simulation Time: 37400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 375, Iterations: 3
Simulation Time: 37500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 376, Iterations: 3
Simulation Time: 37600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 377, Iterations: 3
Simulation Time: 37700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 378, Iterations: 3
Simulation Time: 37800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 379, Iterations: 3
Simulation Time: 37900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 380, Iterations: 3
Simulation Time: 38000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 381, Iterations: 3
Simulation Time: 38100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 382, Iterations: 3
Simulation Time: 38200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 383, Iterations: 3
Simulation Time: 38300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 384, Iterations: 3
Simulation Time: 38400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 385, Iterations: 3
Simulation Time: 38500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 386, Iterations: 3
Simulation Time: 38600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 387, Iterations: 3
Simulation Time: 38700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 388, Iterations: 3
Simulation Time: 38800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 389, Iterations: 3
Simulation Time: 38900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 390, Iterations: 3
Simulation Time: 39000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 391, Iterations: 3
Simulation Time: 39100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 392, Iterations: 3
Simulation Time: 39200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 393, Iterations: 3
Simulation Time: 39300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 394, Iterations: 3
Simulation Time: 39400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 395, Iterations: 3
Simulation Time: 39500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 396, Iterations: 3
Simulation Time: 39600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 397, Iterations: 3
Simulation Time: 39700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 398, Iterations: 3
Simulation Time: 39800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 399, Iterations: 3
Simulation Time: 39900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 400, Iterations: 3
Simulation Time: 40000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 401, Iterations: 3
Simulation Time: 40100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 402, Iterations: 3
Simulation Time: 40200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 403, Iterations: 3
Simulation Time: 40300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 404, Iterations: 3
Simulation Time: 40400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 405, Iterations: 3
Simulation Time: 40500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 406, Iterations: 3
Simulation Time: 40600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 407, Iterations: 3
Simulation Time: 40700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 408, Iterations: 3
Simulation Time: 40800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 409, Iterations: 3
Simulation Time: 40900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 410, Iterations: 3
Simulation Time: 41000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 411, Iterations: 3
Simulation Time: 41100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 412, Iterations: 3
Simulation Time: 41200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 413, Iterations: 3
Simulation Time: 41300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 414, Iterations: 3
Simulation Time: 41400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 415, Iterations: 3
Simulation Time: 41500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 416, Iterations: 3
Simulation Time: 41600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 417, Iterations: 3
Simulation Time: 41700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 418, Iterations: 3
Simulation Time: 41800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 419, Iterations: 3
Simulation Time: 41900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 420, Iterations: 3
Simulation Time: 42000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 421, Iterations: 3
Simulation Time: 42100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 422, Iterations: 3
Simulation Time: 42200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 423, Iterations: 3
Simulation Time: 42300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 424, Iterations: 3
Simulation Time: 42400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 425, Iterations: 3
Simulation Time: 42500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 426, Iterations: 3
Simulation Time: 42600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 427, Iterations: 3
Simulation Time: 42700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 428, Iterations: 3
Simulation Time: 42800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 429, Iterations: 3
Simulation Time: 42900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 430, Iterations: 3
Simulation Time: 43000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 431, Iterations: 3
Simulation Time: 43100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 432, Iterations: 3
Simulation Time: 43200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 433, Iterations: 3
Simulation Time: 43300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 434, Iterations: 3
Simulation Time: 43400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 435, Iterations: 3
Simulation Time: 43500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 436, Iterations: 3
Simulation Time: 43600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 437, Iterations: 3
Simulation Time: 43700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 438, Iterations: 3
Simulation Time: 43800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 439, Iterations: 3
Simulation Time: 43900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 440, Iterations: 3
Simulation Time: 44000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 441, Iterations: 3
Simulation Time: 44100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 442, Iterations: 3
Simulation Time: 44200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 443, Iterations: 3
Simulation Time: 44300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 444, Iterations: 3
Simulation Time: 44400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 445, Iterations: 3
Simulation Time: 44500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 446, Iterations: 3
Simulation Time: 44600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 447, Iterations: 3
Simulation Time: 44700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 448, Iterations: 3
Simulation Time: 44800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 449, Iterations: 3
Simulation Time: 44900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 450, Iterations: 3
Simulation Time: 45000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 451, Iterations: 3
Simulation Time: 45100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 452, Iterations: 3
Simulation Time: 45200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 453, Iterations: 3
Simulation Time: 45300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 454, Iterations: 3
Simulation Time: 45400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 455, Iterations: 3
Simulation Time: 45500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 456, Iterations: 3
Simulation Time: 45600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 457, Iterations: 3
Simulation Time: 45700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 458, Iterations: 3
Simulation Time: 45800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 459, Iterations: 3
Simulation Time: 45900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 460, Iterations: 3
Simulation Time: 46000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 461, Iterations: 3
Simulation Time: 46100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 462, Iterations: 3
Simulation Time: 46200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 463, Iterations: 3
Simulation Time: 46300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 464, Iterations: 3
Simulation Time: 46400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 465, Iterations: 3
Simulation Time: 46500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 466, Iterations: 3
Simulation Time: 46600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 467, Iterations: 3
Simulation Time: 46700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 468, Iterations: 3
Simulation Time: 46800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 469, Iterations: 3
Simulation Time: 46900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 470, Iterations: 3
Simulation Time: 47000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 471, Iterations: 3
Simulation Time: 47100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 472, Iterations: 3
Simulation Time: 47200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 473, Iterations: 3
Simulation Time: 47300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 474, Iterations: 3
Simulation Time: 47400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 475, Iterations: 3
Simulation Time: 47500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 476, Iterations: 3
Simulation Time: 47600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 477, Iterations: 3
Simulation Time: 47700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 478, Iterations: 3
Simulation Time: 47800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 479, Iterations: 3
Simulation Time: 47900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 480, Iterations: 3
Simulation Time: 48000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 481, Iterations: 3
Simulation Time: 48100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 482, Iterations: 3
Simulation Time: 48200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 483, Iterations: 3
Simulation Time: 48300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 484, Iterations: 3
Simulation Time: 48400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 485, Iterations: 3
Simulation Time: 48500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 486, Iterations: 3
Simulation Time: 48600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 487, Iterations: 3
Simulation Time: 48700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 488, Iterations: 3
Simulation Time: 48800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 489, Iterations: 3
Simulation Time: 48900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 490, Iterations: 3
Simulation Time: 49000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 491, Iterations: 3
Simulation Time: 49100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 492, Iterations: 3
Simulation Time: 49200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 493, Iterations: 3
Simulation Time: 49300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 494, Iterations: 3
Simulation Time: 49400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 495, Iterations: 3
Simulation Time: 49500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 496, Iterations: 3
Simulation Time: 49600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 497, Iterations: 3
Simulation Time: 49700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 498, Iterations: 3
Simulation Time: 49800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 499, Iterations: 3
Simulation Time: 49900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 500, Iterations: 3
Simulation Time: 50000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 501, Iterations: 3
Simulation Time: 50100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 502, Iterations: 3
Simulation Time: 50200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 503, Iterations: 3
Simulation Time: 50300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 504, Iterations: 3
Simulation Time: 50400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 505, Iterations: 3
Simulation Time: 50500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 506, Iterations: 3
Simulation Time: 50600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 507, Iterations: 3
Simulation Time: 50700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 508, Iterations: 3
Simulation Time: 50800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 509, Iterations: 3
Simulation Time: 50900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 510, Iterations: 3
Simulation Time: 51000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 511, Iterations: 4
Simulation Time: 51100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 512, Iterations: 4
Simulation Time: 51200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 513, Iterations: 4
Simulation Time: 51300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 514, Iterations: 4
Simulation Time: 51400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 515, Iterations: 4
Simulation Time: 51500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 516, Iterations: 4
Simulation Time: 51600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 517, Iterations: 4
Simulation Time: 51700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 518, Iterations: 4
Simulation Time: 51800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 519, Iterations: 4
Simulation Time: 51900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 520, Iterations: 4
Simulation Time: 52000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 521, Iterations: 4
Simulation Time: 52100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 522, Iterations: 4
Simulation Time: 52200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 523, Iterations: 4
Simulation Time: 52300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 524, Iterations: 4
Simulation Time: 52400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 525, Iterations: 4
Simulation Time: 52500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 526, Iterations: 4
Simulation Time: 52600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 527, Iterations: 4
Simulation Time: 52700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 528, Iterations: 4
Simulation Time: 52800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 529, Iterations: 4
Simulation Time: 52900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 530, Iterations: 4
Simulation Time: 53000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 531, Iterations: 4
Simulation Time: 53100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 532, Iterations: 4
Simulation Time: 53200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 533, Iterations: 4
Simulation Time: 53300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 534, Iterations: 4
Simulation Time: 53400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 535, Iterations: 4
Simulation Time: 53500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 536, Iterations: 4
Simulation Time: 53600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 537, Iterations: 4
Simulation Time: 53700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 538, Iterations: 4
Simulation Time: 53800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 539, Iterations: 4
Simulation Time: 53900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 540, Iterations: 4
Simulation Time: 54000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 541, Iterations: 4
Simulation Time: 54100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 542, Iterations: 4
Simulation Time: 54200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 543, Iterations: 4
Simulation Time: 54300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 544, Iterations: 4
Simulation Time: 54400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 545, Iterations: 4
Simulation Time: 54500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 546, Iterations: 4
Simulation Time: 54600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 547, Iterations: 4
Simulation Time: 54700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 548, Iterations: 4
Simulation Time: 54800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 549, Iterations: 4
Simulation Time: 54900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 550, Iterations: 4
Simulation Time: 55000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 551, Iterations: 4
Simulation Time: 55100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 552, Iterations: 4
Simulation Time: 55200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 553, Iterations: 4
Simulation Time: 55300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 554, Iterations: 4
Simulation Time: 55400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 555, Iterations: 4
Simulation Time: 55500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 556, Iterations: 4
Simulation Time: 55600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 557, Iterations: 4
Simulation Time: 55700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 558, Iterations: 4
Simulation Time: 55800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 559, Iterations: 4
Simulation Time: 55900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 560, Iterations: 4
Simulation Time: 56000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 561, Iterations: 4
Simulation Time: 56100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 562, Iterations: 4
Simulation Time: 56200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 563, Iterations: 4
Simulation Time: 56300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 564, Iterations: 4
Simulation Time: 56400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 565, Iterations: 4
Simulation Time: 56500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 566, Iterations: 4
Simulation Time: 56600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 567, Iterations: 4
Simulation Time: 56700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 568, Iterations: 4
Simulation Time: 56800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 569, Iterations: 4
Simulation Time: 56900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 570, Iterations: 4
Simulation Time: 57000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 571, Iterations: 4
Simulation Time: 57100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 572, Iterations: 4
Simulation Time: 57200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 573, Iterations: 4
Simulation Time: 57300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 574, Iterations: 4
Simulation Time: 57400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 575, Iterations: 4
Simulation Time: 57500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 576, Iterations: 4
Simulation Time: 57600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 577, Iterations: 3
Simulation Time: 57700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 578, Iterations: 3
Simulation Time: 57800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 579, Iterations: 3
Simulation Time: 57900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 580, Iterations: 3
Simulation Time: 58000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 581, Iterations: 3
Simulation Time: 58100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 582, Iterations: 3
Simulation Time: 58200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 583, Iterations: 3
Simulation Time: 58300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 584, Iterations: 3
Simulation Time: 58400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 585, Iterations: 3
Simulation Time: 58500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 586, Iterations: 3
Simulation Time: 58600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 587, Iterations: 3
Simulation Time: 58700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 588, Iterations: 3
Simulation Time: 58800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 589, Iterations: 3
Simulation Time: 58900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 590, Iterations: 3
Simulation Time: 59000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 591, Iterations: 3
Simulation Time: 59100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 592, Iterations: 3
Simulation Time: 59200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 593, Iterations: 3
Simulation Time: 59300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 594, Iterations: 3
Simulation Time: 59400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 595, Iterations: 3
Simulation Time: 59500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 596, Iterations: 3
Simulation Time: 59600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 597, Iterations: 3
Simulation Time: 59700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 598, Iterations: 3
Simulation Time: 59800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 599, Iterations: 3
Simulation Time: 59900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 600, Iterations: 3
Simulation Time: 60000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 601, Iterations: 3
Simulation Time: 60100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 602, Iterations: 3
Simulation Time: 60200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 603, Iterations: 3
Simulation Time: 60300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 604, Iterations: 3
Simulation Time: 60400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 605, Iterations: 3
Simulation Time: 60500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 606, Iterations: 3
Simulation Time: 60600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 607, Iterations: 3
Simulation Time: 60700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 608, Iterations: 3
Simulation Time: 60800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 609, Iterations: 3
Simulation Time: 60900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 610, Iterations: 3
Simulation Time: 61000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 611, Iterations: 3
Simulation Time: 61100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 612, Iterations: 3
Simulation Time: 61200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 613, Iterations: 3
Simulation Time: 61300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 614, Iterations: 3
Simulation Time: 61400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 615, Iterations: 3
Simulation Time: 61500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 616, Iterations: 3
Simulation Time: 61600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 617, Iterations: 3
Simulation Time: 61700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 618, Iterations: 3
Simulation Time: 61800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 619, Iterations: 3
Simulation Time: 61900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 620, Iterations: 3
Simulation Time: 62000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 621, Iterations: 3
Simulation Time: 62100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 622, Iterations: 3
Simulation Time: 62200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 623, Iterations: 3
Simulation Time: 62300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 624, Iterations: 3
Simulation Time: 62400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 625, Iterations: 3
Simulation Time: 62500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 626, Iterations: 3
Simulation Time: 62600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 627, Iterations: 3
Simulation Time: 62700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 628, Iterations: 3
Simulation Time: 62800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 629, Iterations: 3
Simulation Time: 62900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 630, Iterations: 3
Simulation Time: 63000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 631, Iterations: 3
Simulation Time: 63100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 632, Iterations: 3
Simulation Time: 63200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 633, Iterations: 3
Simulation Time: 63300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 634, Iterations: 3
Simulation Time: 63400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 635, Iterations: 3
Simulation Time: 63500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 636, Iterations: 3
Simulation Time: 63600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 637, Iterations: 3
Simulation Time: 63700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 638, Iterations: 3
Simulation Time: 63800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 639, Iterations: 3
Simulation Time: 63900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 640, Iterations: 3
Simulation Time: 64000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 641, Iterations: 3
Simulation Time: 64100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 642, Iterations: 3
Simulation Time: 64200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 643, Iterations: 3
Simulation Time: 64300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 644, Iterations: 3
Simulation Time: 64400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 645, Iterations: 3
Simulation Time: 64500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 646, Iterations: 3
Simulation Time: 64600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 647, Iterations: 3
Simulation Time: 64700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 648, Iterations: 3
Simulation Time: 64800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 649, Iterations: 3
Simulation Time: 64900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 650, Iterations: 3
Simulation Time: 65000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 651, Iterations: 3
Simulation Time: 65100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 652, Iterations: 3
Simulation Time: 65200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 653, Iterations: 3
Simulation Time: 65300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 654, Iterations: 3
Simulation Time: 65400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 655, Iterations: 3
Simulation Time: 65500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 656, Iterations: 3
Simulation Time: 65600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 657, Iterations: 3
Simulation Time: 65700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 658, Iterations: 3
Simulation Time: 65800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 659, Iterations: 3
Simulation Time: 65900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 660, Iterations: 3
Simulation Time: 66000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 661, Iterations: 3
Simulation Time: 66100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 662, Iterations: 3
Simulation Time: 66200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 663, Iterations: 3
Simulation Time: 66300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 664, Iterations: 3
Simulation Time: 66400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 665, Iterations: 3
Simulation Time: 66500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 666, Iterations: 3
Simulation Time: 66600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 667, Iterations: 3
Simulation Time: 66700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 668, Iterations: 3
Simulation Time: 66800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 669, Iterations: 3
Simulation Time: 66900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 670, Iterations: 3
Simulation Time: 67000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 671, Iterations: 3
Simulation Time: 67100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 672, Iterations: 3
Simulation Time: 67200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 673, Iterations: 3
Simulation Time: 67300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 674, Iterations: 3
Simulation Time: 67400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 675, Iterations: 3
Simulation Time: 67500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 676, Iterations: 3
Simulation Time: 67600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 677, Iterations: 3
Simulation Time: 67700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 678, Iterations: 3
Simulation Time: 67800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 679, Iterations: 3
Simulation Time: 67900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 680, Iterations: 3
Simulation Time: 68000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 681, Iterations: 3
Simulation Time: 68100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 682, Iterations: 3
Simulation Time: 68200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 683, Iterations: 3
Simulation Time: 68300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 684, Iterations: 3
Simulation Time: 68400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 685, Iterations: 3
Simulation Time: 68500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 686, Iterations: 3
Simulation Time: 68600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 687, Iterations: 3
Simulation Time: 68700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 688, Iterations: 3
Simulation Time: 68800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 689, Iterations: 3
Simulation Time: 68900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 690, Iterations: 3
Simulation Time: 69000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 691, Iterations: 3
Simulation Time: 69100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 692, Iterations: 3
Simulation Time: 69200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 693, Iterations: 3
Simulation Time: 69300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 694, Iterations: 3
Simulation Time: 69400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 695, Iterations: 3
Simulation Time: 69500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 696, Iterations: 3
Simulation Time: 69600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 697, Iterations: 3
Simulation Time: 69700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 698, Iterations: 3
Simulation Time: 69800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 699, Iterations: 3
Simulation Time: 69900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 700, Iterations: 3
Simulation Time: 70000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 701, Iterations: 3
Simulation Time: 70100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 702, Iterations: 3
Simulation Time: 70200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 703, Iterations: 3
Simulation Time: 70300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 704, Iterations: 3
Simulation Time: 70400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 705, Iterations: 3
Simulation Time: 70500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 706, Iterations: 3
Simulation Time: 70600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 707, Iterations: 3
Simulation Time: 70700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 708, Iterations: 3
Simulation Time: 70800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 709, Iterations: 3
Simulation Time: 70900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 710, Iterations: 3
Simulation Time: 71000.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 711, Iterations: 3
Simulation Time: 71100.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 712, Iterations: 3
Simulation Time: 71200.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 713, Iterations: 3
Simulation Time: 71300.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 714, Iterations: 3
Simulation Time: 71400.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 715, Iterations: 3
Simulation Time: 71500.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 716, Iterations: 3
Simulation Time: 71600.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 717, Iterations: 3
Simulation Time: 71700.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 718, Iterations: 3
Simulation Time: 71800.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 719, Iterations: 3
Simulation Time: 71900.0 s of 72000 s
Step: Swell-deswell-reswell | Increment: 720, Iterations: 3
Simulation Time: 72000.0 s of 72000 s
-----------------------------------------
End computation
------------------------------------------
Elapsed real time: 0:00:47.732805
------------------------------------------
# Set up font size, initialize colors array
font = {'size' : 14}
plt.rc('font', **font)
#
prop_cycle = plt.rcParams['axes.prop_cycle']
colors = prop_cycle.by_key()['color']
# Create figure for polymer volume fraction phi = (1+c)**(-1) versus normalized position x/L0 at different times:
#
fig = plt.figure()
ax=fig.gca()
#---------------------------------------------------------------------------------------------
for i in range(len(profileTimes)):
ax.plot(xPath/2.5, (1+C_res[:, i])**(-1.0), label=r"$t =$ {} h".format(profileTimes[i]/3600),\
linewidth=2.0)
#---------------------------------------------------------------------------------------------
ax.set_xlabel(r"$x$-coordinate along $y=0$", size=14)
ax.set_ylabel(r"Polymer volume fraction $\phi$", size=14)
ax.set_ylim(0.15, 0.45)
ax.set_xlim(0.0,1.0)
#
plt.grid(linestyle="--", linewidth=0.5, color='b')
#
from matplotlib.ticker import AutoMinorLocator,FormatStrFormatter
ax.xaxis.set_minor_locator(AutoMinorLocator())
ax.yaxis.set_minor_locator(AutoMinorLocator())
#
plt.legend(bbox_to_anchor=(1.0, 1.0))
import matplotlib.ticker as ticker
# ax.xaxis.set_major_formatter(ticker.FormatStrFormatter('%0.0f'))
# Save figure
fig = plt.gcf()
fig.set_size_inches(8,5)
plt.tight_layout()
plt.savefig("results/gel_pe_phi_time_profiles.png", dpi=600)
