Commit a8b9a7bb by Jigyasa Watwani

fully working 1-D problem with remeshing and extend functionality

parent 3a76deec
...@@ -15,6 +15,6 @@ ...@@ -15,6 +15,6 @@
"noise_level": 0.0, "noise_level": 0.0,
"timestep": 0.01, "timestep": 0.01,
"savetime": 0.1, "savetime": 0.1,
"maxtime": 100 "maxtime": 10
} }
\ No newline at end of file
...@@ -36,17 +36,28 @@ else: ...@@ -36,17 +36,28 @@ else:
savesteps = round(parameters['maxtime']/parameters['savetime']) savesteps = round(parameters['maxtime']/parameters['savetime'])
mesh = df.IntervalMesh(parameters['resolution'], # read geometry and topology at the last timepoint from h5 file
- parameters['system_size']/2,
parameters['system_size']/2)
# read geometry from h5 file
var = 'density' var = 'density'
h5 = h5py.File( '%s_%s.h5' % (parameters['timestamp'], var,), 'r+') h5 = h5py.File( '%s_%s.h5' % (parameters['timestamp'], var,), 'r+')
geometry = np.array(h5['%s/%s_%d/mesh/geometry'%(var,var, savesteps)]) geometry = np.array(h5['%s/%s_%d/mesh/geometry'%(var,var,savesteps)])[:, 0]
mesh.coordinates()[:,0] = geometry[:,0] topology = np.array(h5['%s/%s_%d/mesh/topology'%(var,var,savesteps)])
# create mesh with this geometry and topology
mesh = df.Mesh()
editor = df.MeshEditor()
editor.open(mesh, 'interval', 1, 1)
# initialise vertices and cells
editor.init_vertices(len(geometry))
editor.init_cells(len(topology))
# add vertices
for j in range(len(geometry)):
editor.add_vertex(j, [geometry[j]])
# add cells
for j in range(len(topology)):
editor.add_cell(j, topology[j])
editor.close()
# Read data # Read field values at the last time point
SFS = df.FunctionSpace(mesh, 'P', 1) SFS = df.FunctionSpace(mesh, 'P', 1)
VFS = df.VectorFunctionSpace(mesh, 'P', 1) VFS = df.VectorFunctionSpace(mesh, 'P', 1)
initu, initv, initrho = df.Function(VFS), df.Function(VFS), df.Function(SFS) initu, initv, initrho = df.Function(VFS), df.Function(VFS), df.Function(SFS)
...@@ -71,5 +82,4 @@ if extend: ...@@ -71,5 +82,4 @@ if extend:
with open(parametersFile, "w") as jsonFile: with open(parametersFile, "w") as jsonFile:
json.dump(parameters, jsonFile, indent=4, sort_keys=True) json.dump(parameters, jsonFile, indent=4, sort_keys=True)
from viz_tissue import visualize
visualize(parameters, DIR='') visualize(parameters, DIR='')
\ No newline at end of file
...@@ -13,22 +13,20 @@ class Tissue(object): ...@@ -13,22 +13,20 @@ class Tissue(object):
# if no mesh is given as initial condition, create one # if no mesh is given as initial condition, create one
if mesh is None: if mesh is None:
self.mesh = df.IntervalMesh(self.resolution, self.mesh = df.IntervalMesh(self.resolution, -self.system_size/2, self.system_size/2)
-self.system_size/2,
self.system_size/2)
else: else:
self.mesh = mesh self.mesh = mesh
def setup(self):
scalar_element = df.FiniteElement('P', self.mesh.ufl_cell(), 1) scalar_element = df.FiniteElement('P', self.mesh.ufl_cell(), 1)
vector_element = df.VectorElement('P', self.mesh.ufl_cell(), 1) vector_element = df.VectorElement('P', self.mesh.ufl_cell(), 1)
# u, v, rho # u, v, rho
mixed_element = df.MixedElement([vector_element, self.mixed_element = df.MixedElement([vector_element, vector_element, scalar_element])
vector_element,
scalar_element])
# define function space with this mixed element # define function space with this mixed element
self.function_space = df.FunctionSpace(self.mesh, mixed_element) self.function_space = df.FunctionSpace(self.mesh, self.mixed_element)
# boundary condition on rho # boundary condition on rho
if self.zero_rho_boundary: if self.zero_rho_boundary:
...@@ -44,7 +42,6 @@ class Tissue(object): ...@@ -44,7 +42,6 @@ class Tissue(object):
self.function0 = df.Function(self.function_space) self.function0 = df.Function(self.function_space)
self.function = df.Function(self.function_space) self.function = df.Function(self.function_space)
def advection(self, conc, vel, tconc): def advection(self, conc, vel, tconc):
return df.inner(df.div(vel * conc), tconc) return df.inner(df.div(vel * conc), tconc)
...@@ -56,10 +53,18 @@ class Tissue(object): ...@@ -56,10 +53,18 @@ class Tissue(object):
return (-self.lamda * rho /(rho + self.saturation_rho) return (-self.lamda * rho /(rho + self.saturation_rho)
* df.Identity(1)) * df.Identity(1))
def epsilon(self, v): def epsilon(self, v):
return df.sym(df.nabla_grad(v)) return df.sym(df.nabla_grad(v))
def refine_mesh(self, mesh):
cell_markers = df.MeshFunction("bool", mesh, mesh.topology().dim())
cell_markers.set_all(False)
for cell in df.cells(mesh):
if cell.h() > 0.1: # if the width of the cell> 0.1, refine
cell_markers[cell] = True
mesh = df.refine(mesh, cell_markers)
return mesh
def passive_stress(self, u, v): def passive_stress(self, u, v):
eps_u, eps_v = self.epsilon(u), self.epsilon(v) eps_u, eps_v = self.epsilon(u), self.epsilon(v)
elastic_stress = self.elasticity * eps_u elastic_stress = self.elasticity * eps_u
...@@ -77,10 +82,11 @@ class Tissue(object): ...@@ -77,10 +82,11 @@ class Tissue(object):
) )
def setup_initial_conditions(self): def setup_initial_conditions(self):
# ic for u
zero_vector = df.Constant((0.0,)) zero_vector = df.Constant((0.0,))
u0 = df.interpolate(zero_vector, u0 = df.interpolate(zero_vector,
self.function_space.sub(0).collapse()) self.function_space.sub(0).collapse())
# ic for rho
if self.zero_rho_boundary: if self.zero_rho_boundary:
base_rho = 0.0 base_rho = 0.0
else: else:
...@@ -123,61 +129,81 @@ class Tissue(object): ...@@ -123,61 +129,81 @@ class Tissue(object):
self.form = (uform + vform + rhoform) * df.dx self.form = (uform + vform + rhoform) * df.dx
def solve(self, initu=None, initv=None, initrho=None, initTime=0, extend=False): def solve(self, initu=None, initv=None, initrho=None, initTime=0, extend=False):
# create function space, functions, bc
self.setup()
# create files if they don't exist, open them if they do exist # create files if they don't exist, open them if they do exist
self.uFile = df.XDMFFile( self.uFile = df.XDMFFile('%s_displacement.xdmf' % self.timestamp)
'%s_displacement.xdmf' % self.timestamp) self.vFile = df.XDMFFile('%s_velocity.xdmf' % self.timestamp)
self.vFile = df.XDMFFile( self.rhoFile = df.XDMFFile('%s_density.xdmf' % self.timestamp)
'%s_velocity.xdmf' % self.timestamp)
self.rhoFile = df.XDMFFile(
'%s_density.xdmf' % self.timestamp)
# time-variables # time-variables
self.time = initTime self.time = initTime
savesteps = round(self.savetime/self.timestep) savesteps = round(self.savetime/self.timestep)
maxsteps = round(self.maxtime/self.timestep) maxsteps = round(self.maxtime/self.timestep)
# if not extend, find the initial velocity, save the ic, move the mesh with this initial velocity
if not extend: if not extend:
# find the initial velocity
self.setup_initial_conditions() self.setup_initial_conditions()
# save the ic
u0, v0, rho0 = self.function0.split(deepcopy=True) u0, v0, rho0 = self.function0.split(deepcopy=True)
self.uFile.write_checkpoint(u0, 'displacement', self.time) self.uFile.write_checkpoint(u0, 'displacement', self.time)
self.vFile.write_checkpoint(v0, 'velocity', self.time) self.vFile.write_checkpoint(v0, 'velocity', self.time)
self.rhoFile.write_checkpoint(rho0, 'density', self.time) self.rhoFile.write_checkpoint(rho0, 'density', self.time)
# move the mesh with this initial velocity
dr0 = df.project(v0 * self.timestep, self.function_space.sub(0).collapse()) dr0 = df.project(v0 * self.timestep, self.function_space.sub(0).collapse())
df.ALE.move(self.mesh, dr0) df.ALE.move(self.mesh, dr0)
# if extend, move the mesh with this initial velocity
else: else:
# assign fields to func0
u0 = df.project(initu, self.function_space.sub(0).collapse()) u0 = df.project(initu, self.function_space.sub(0).collapse())
v0 = df.project(initv, self.function_space.sub(1).collapse()) v0 = df.project(initv, self.function_space.sub(1).collapse())
rho0 = df.project(initrho, self.function_space.sub(2).collapse()) rho0 = df.project(initrho, self.function_space.sub(2).collapse())
df.assign(self.function0, [u0, v0, rho0])
# move the mesh with initial v read from file
dr0 = df.project(v0 * self.timestep, self.function_space.sub(0).collapse()) dr0 = df.project(v0 * self.timestep, self.function_space.sub(0).collapse())
df.ALE.move(self.mesh, dr0) df.ALE.move(self.mesh, dr0)
df.assign(self.function0, [u0, v0, rho0])
self.setup_weak_forms() self.setup_weak_forms()
for steps in progressbar.progressbar(range(1, maxsteps + 1)): for steps in progressbar.progressbar(range(1, maxsteps + 1)):
self.time += self.timestep self.time += self.timestep
# solve to get fields at next timestep
df.solve(self.form == 0, self.function, self.bc) df.solve(self.form == 0, self.function, self.bc)
self.function0.assign(self.function)
u, v, rho = self.function0.split(deepcopy=True) # save the fields
u, v, rho = self.function.split(deepcopy=True)
if steps % savesteps == 0: if steps % savesteps == 0:
self.uFile.write_checkpoint(u, 'displacement', self.uFile.write_checkpoint(u, 'displacement', self.time, append=True)
self.time, append=True) self.vFile.write_checkpoint(v, 'velocity', self.time, append=True)
self.vFile.write_checkpoint(v, 'velocity', self.rhoFile.write_checkpoint(rho, 'density', self.time, append=True)
self.time, append=True)
self.rhoFile.write_checkpoint(rho, 'density', # assign the calculated function to a newly-defined function, to be used later
self.time, append=True) old_function = df.Function(self.function_space)
old_function.assign(self.function)
# move mesh # move the mesh
dr = df.project(v * self.timestep, self.function_space.sub(0).collapse()) dr = df.project(v * self.timestep, self.function_space.sub(0).collapse())
df.ALE.move(self.mesh, dr) df.ALE.move(self.mesh, dr)
# refine the mesh
self.mesh = self.refine_mesh(self.mesh)
# new function space, bc, functions
self.setup()
# new function0 should be the old function in the old function space projected on the new function space
self.function0 = df.project(old_function, self.function_space)
# new weak form
self.setup_weak_forms()
self.uFile.close() self.uFile.close()
self.vFile.close() self.vFile.close()
......
Markdown is supported
0% or
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or sign in to comment