217 lines
6.7 KiB
Python
217 lines
6.7 KiB
Python
from dolfin import *
|
|
import matplotlib.pyplot as plt
|
|
import numpy as np
|
|
from common import inout
|
|
from mpi4py import MPI
|
|
import sys
|
|
import os
|
|
#
|
|
# NAVIER STOKES PROBLEM IN THE AORTA with a MONOLITHIC SOLVER
|
|
# THIS SCRIPT INCLUDE THE 0-WINDKESSEL BOUNDARY CONDITION
|
|
#
|
|
# Written by Jeremias Garay L: j.e.garay.labra@rug.nl
|
|
#
|
|
|
|
parameters["std_out_all_processes"] = False
|
|
|
|
def solv_NavierStokes(options):
|
|
|
|
# Assign physical parameters
|
|
rho = Constant(options['density'])
|
|
mu = Constant(options['dynamic_viscosity'])
|
|
otheta = Constant(1-options['theta'])
|
|
theta = Constant(options['theta'])
|
|
Tf = options['Tf']
|
|
dt = options['dt']
|
|
dt_w = options['dt_write']
|
|
|
|
# CREATING THE FILES
|
|
xdmf_u = XDMFFile(options['savepath']+'u.xdmf')
|
|
xdmf_p = XDMFFile(options['savepath']+'p.xdmf')
|
|
# LOADING THE MESH
|
|
mesh = Mesh()
|
|
hdf = HDF5File(mesh.mpi_comm(), options['mesh_path'] , 'r')
|
|
hdf.read(mesh, '/mesh', False)
|
|
bnds = MeshFunction('size_t', mesh , mesh.topology().dim() - 1)
|
|
hdf.read(bnds, '/boundaries')
|
|
|
|
# DEFINE FUNCTION SPACES
|
|
if options['fem_space'] == 'p2p1':
|
|
V = VectorElement('Lagrange', mesh.ufl_cell(), 2)
|
|
Q = FiniteElement('Lagrange', mesh.ufl_cell(), 1)
|
|
pspg = False
|
|
elif options['fem_space'] == 'p1p1':
|
|
V = VectorElement('Lagrange', mesh.ufl_cell(), 1)
|
|
Q = FiniteElement('Lagrange', mesh.ufl_cell(), 1)
|
|
pspg = True
|
|
TH = V * Q
|
|
W = FunctionSpace(mesh, TH)
|
|
n = FacetNormal(mesh)
|
|
ds = Measure("ds", subdomain_data=bnds)
|
|
v, q = TestFunctions(W)
|
|
# Boundary Conditions
|
|
BCS = []
|
|
bc = options['boundary_conditions']
|
|
# For Windkessel implementation
|
|
flows = {}
|
|
pii0 = {}
|
|
pii = {}
|
|
press = {}
|
|
alpha = {}
|
|
beta = {}
|
|
gamma = {}
|
|
Windkvar = {}
|
|
windkessel = False
|
|
|
|
for nbc in range(len(bc)):
|
|
nid = bc[nbc]['id']
|
|
if bc[nbc]['type'] == 'dirichlet':
|
|
if rank==0:
|
|
print('Adding Dirichlet BC at boundary ', nid)
|
|
|
|
val = bc[nbc]['value']
|
|
if isinstance(val, (int, float)):
|
|
val = Constant(val)
|
|
BCS.append(DirichletBC(W.sub(0), val, bnds, nid))
|
|
else:
|
|
params = bc[nbc]['parameters'] if 'parameters' in bc[nbc] else dict()
|
|
inflow = Expression(val, degree=3, **params)
|
|
BCS.append(DirichletBC(W.sub(0), inflow, bnds, nid))
|
|
|
|
if bc[nbc]['type'] == 'windkessel':
|
|
windkessel = True
|
|
if rank==0:
|
|
print('Adding Windkessel BC at boundary ',nid)
|
|
[R_p,C,R_d] = bc[nbc]['value']
|
|
# coeficients
|
|
alpha[nid] = R_d*C/(R_d*C + dt)
|
|
beta[nid] = R_d*(1-alpha[nid])
|
|
gamma[nid] = R_p + beta[nid]
|
|
# dynamical terms
|
|
if C ==0:
|
|
flows[nid] = Constant(0)
|
|
press[nid] = Constant(R_p*0)
|
|
else:
|
|
flows[nid] = Constant(0)
|
|
pii0[nid] = Constant(bc[nbc]['p0'][0]*bc[nbc]['p0'][1])
|
|
pii[nid] = Constant(alpha[nid]*pii0[nid] + beta[nid]*flows[nid])
|
|
press[nid] = Constant(gamma[nid]*flows[nid] + alpha[nid]*pii0[nid])
|
|
|
|
Windkvar[nid] = dt*press[nid]*dot(v,n)*ds(nid)
|
|
|
|
u, p = TrialFunctions(W)
|
|
w = Function(W)
|
|
h = CellDiameter(mesh)
|
|
|
|
u0, p0 = w.split()
|
|
u_ = theta*u + otheta*u0
|
|
p_ = theta*p + otheta*p0
|
|
# The variational formulation
|
|
# Mass matrix
|
|
F = (
|
|
(rho/dt)*dot(u - u0, v)*dx
|
|
+ mu*inner(grad(u_), grad(v))*dx
|
|
- p_*div(v)*dx + q*div(u)*dx
|
|
+ rho*dot(grad(u_)*u0, v)*dx
|
|
)
|
|
|
|
# Stabilization Terms
|
|
if options['stabilization']['temam']:
|
|
if rank==0:
|
|
print('Addint Temam stabilization term')
|
|
F += 0.5*rho*div(u0)*dot(u_, v)*dx
|
|
|
|
if pspg:
|
|
if rank==0:
|
|
print('Adding PSPG stabilization term')
|
|
eps = Constant(options['stabilization']['eps'])
|
|
F += eps/mu*h**2*inner(grad(p_), grad(q))*dx
|
|
|
|
if options['stabilization']['forced_normal']['enabled']:
|
|
gparam = options['stabilization']['forced_normal']['gamma']
|
|
for nid in options['stabilization']['forced_normal']['boundaries']:
|
|
if rank==0:
|
|
print('Forcing normal velocity in border ', nid)
|
|
ut = u - n*dot(u,n)
|
|
vt = v - n*dot(v,n)
|
|
F += gparam*dot(ut,vt)*ds(nid)
|
|
|
|
|
|
if len(options['stabilization']['backflow_boundaries'])>0:
|
|
def abs_n(x):
|
|
return 0.5*(x - abs(x))
|
|
for nk in options['stabilization']['backflow_boundaries']:
|
|
if rank==0:
|
|
print('adding backflow stabilization in border number:' + str(nk))
|
|
F -= 0.5*rho*abs_n(dot(u0, n))*dot(u_, v)*ds(nk)
|
|
|
|
a = lhs(F)
|
|
L = rhs(F)
|
|
|
|
if windkessel:
|
|
for nid in Windkvar.keys():
|
|
L = L - Windkvar[nid]
|
|
# The static part of the matrix
|
|
A = assemble(a)
|
|
u, p = w.split()
|
|
u.rename('velocity', 'u')
|
|
p.rename('pressure', 'p')
|
|
ind = 0
|
|
t = dt
|
|
checkcicle = int(options['checkpoint_dt']/options['dt'])
|
|
writecicle = int(options['checkpoint_dt']/options['dt_write'])
|
|
|
|
while t<=Tf+dt:
|
|
if windkessel:
|
|
for k in flows.keys():
|
|
flows[k].assign(assemble(dot(u0,n)*ds(k)))
|
|
pii0[k].assign(pii[k])
|
|
pii[k].assign(alpha[k]*pii0[k] + beta[k]*flows[k])
|
|
press[k].assign(gamma[k]*flows[k] + alpha[k]*pii0[k])
|
|
|
|
# To solve
|
|
assemble(a, tensor=A)
|
|
b = assemble(L)
|
|
[bcs.apply(A, b) for bcs in BCS]
|
|
solve(A, w.vector(), b , 'lu')
|
|
ind += 1
|
|
if options['write_xdmf']:
|
|
if np.mod(ind,writecicle)<0.1 or ind==1:
|
|
xdmf_u.write(u, t)
|
|
xdmf_p.write(p, t)
|
|
|
|
if np.mod(ind,checkcicle)<0.1 or ind==1:
|
|
if options['write_checkpoint']:
|
|
checkpath = options['savepath'] +'checkpoint/{i}/'.format(i=ind)
|
|
comm = u.function_space().mesh().mpi_comm()
|
|
inout.write_HDF5_data(comm, checkpath + '/u.h5', u, '/u', t=t)
|
|
inout.write_HDF5_data(comm, checkpath + '/p.h5', p, '/p', t=t)
|
|
|
|
inflow.t = t
|
|
t += dt
|
|
assign(u0, w.sub(0))
|
|
|
|
|
|
|
|
if __name__ == '__main__':
|
|
|
|
comm = MPI.COMM_WORLD
|
|
size = comm.Get_size()
|
|
rank = comm.Get_rank()
|
|
|
|
if len(sys.argv) > 1:
|
|
if os.path.exists(sys.argv[1]):
|
|
inputfile = sys.argv[1]
|
|
if rank==0:
|
|
print('Found input file ' + inputfile)
|
|
else:
|
|
raise Exception('Command line arg given but input file does not exist:'
|
|
' {}'.format(sys.argv[1]))
|
|
else:
|
|
raise Exception('An input file is required as argument!')
|
|
|
|
|
|
options = inout.read_parameters(inputfile)
|
|
solv_NavierStokes(options)
|
|
|