Newer
Older
#coding=utf8
########################################################################
### ###
### Created by Martin Genet, 2016 ###
### ###
### École Polytechnique, Palaiseau, France ###
### ###
########################################################################
import dolfin
import glob
import numpy
import os
import myFEniCSPythonLibrary as myFEniCS
import myVTKPythonLibrary as myVTK
########################################################################
def print_str(tab, string):
print " | "*tab + string
def print_var(tab, name, val):
print " | "*tab + name + " = " + str(val)
def print_sci(tab, name, val):
print " | "*tab + name.ljust(12) + " = " + format(val,".4e")
def fedic(
working_folder,
working_basename,
images_folder,
images_basename,
mesh_folder,
mesh_basename,

Martin Genet
committed
images_quadrature=None,

Martin Genet
committed
penalty=0.9,
tangent_type="Idef", # Idef, Idef-wHess, Idef_old, Idef_old-wHess, Iref

Martin Genet
committed
continue_after_fail=0):
assert os.path.exists(mesh_filename), "No mesh in "+mesh_filename+". Aborting."
mesh = dolfin.Mesh(mesh_filename)
dX = dolfin.dx(mesh)
print_str(tab,"Computing quadrature degree for images…")
ref_image_filename = images_folder+"/"+images_basename+"_"+str(images_k_ref).zfill(images_zfill)+".vti"

Martin Genet
committed
if (images_quadrature is None):
images_quadrature = myFEniCS.compute_quadrature_degree(
image_filename=ref_image_filename,
mesh=mesh,
verbose=0)
print_var(tab+1,"images_quadrature",images_quadrature)

Martin Genet
committed
image_filename=ref_image_filename,
verbose=0)
assert (images_dimension in (2,3)), "images_dimension must be 2 or 3. Aborting."
cell=mesh.ufl_cell(),
quad_scheme="default")
ve = dolfin.VectorElement(
family="Quadrature",
cell=mesh.ufl_cell(),
quad_scheme="default")

Martin Genet
committed
te = dolfin.TensorElement(
family="Quadrature",
cell=mesh.ufl_cell(),
degree=images_quadrature,
quad_scheme="default")
te._quad_scheme = "default"
for k in xrange(images_dimension**2):
te.sub_elements()[k]._quad_scheme = "default"
Iref_norm = (dolfin.assemble(Iref**2 * dX, form_compiler_parameters={'quadrature_degree':images_quadrature})/mesh_volume)**(1./2)
assert (Iref_norm > 0.), "Iref_norm = "+str(Iref_norm)+" <= 0. Aborting."
vfs = dolfin.VectorFunctionSpace(
mesh=mesh,
family="Lagrange",
degree=mesh_degree)
U_old = dolfin.Function(

Martin Genet
committed
DU_old = dolfin.Function(
vfs,
name="previous displacement increment")
dU_ = dolfin.TrialFunction(vfs)
dV_ = dolfin.TestFunction(vfs)

Martin Genet
committed
DDIdef = myFEniCS.ExprHessDefIm2(
U=U,
filename=ref_image_filename,
element=te)

Martin Genet
committed
DDIdef = myFEniCS.ExprHessDefIm3(
U=U,
filename=ref_image_filename,
element=te)
if (tangent_type == "Iref"):
a = penalty * dolfin.inner(dolfin.dot(DIref, dU_),

Martin Genet
committed
dolfin.dot(DIref, dV_)) * dX
elif (tangent_type.startswith("Idef")):

Martin Genet
committed
dolfin.dot(DIdef, dV_)) * dX
if ("-wHess" in tangent_type):
a += penalty * (Idef-Iref) * dolfin.inner(dolfin.dot(DDIdef, dU_),
dV_) * dX
b = - penalty * (Idef-Iref) * dolfin.dot(DIdef, dV_) * dX

Martin Genet
committed
dolfin.grad(dV_)) * dX
b += - (1.-penalty) * dolfin.inner(dolfin.grad( U ),
dolfin.grad(dV_)) * dX
b0 = dolfin.inner(Iref, dolfin.dot(DIref, dV_)) * dX
B0 = dolfin.assemble(b0, form_compiler_parameters={'quadrature_degree':images_quadrature})
assert (B0_norm > 0.), "B0_norm = "+str(B0_norm)+" <= 0. Aborting."
A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})

Martin Genet
committed
n_frames = len(glob.glob(images_folder+"/"+images_basename+"_"+"?"*images_zfill+".vti"))
#n_frames = 2
assert (n_frames > 1), "n_frames = "+str(n_frames)+" <= 1. Aborting."
assert (abs(images_k_ref) < n_frames), "abs(images_k_ref) = "+str(images_k_ref)+" >= n_frames. Aborting."
if not os.path.exists(working_folder):
os.mkdir(working_folder)
pvd_basename = working_folder+"/"+working_basename+"_"
file_pvd = dolfin.File(pvd_basename+".pvd")
for vtu_filename in glob.glob(pvd_basename+"*.vtu"):
os.remove(vtu_filename)
file_pvd << (U, float(images_k_ref))
for forward_or_backward in ["forward", "backward"]: # not sure if this still works…
print_var(tab,"forward_or_backward",forward_or_backward)
if (forward_or_backward == "backward"):
U.vector()[:] = 0.

Martin Genet
committed
U_old.vector()[:] = 0.
DU_old.vector()[:] = 0.
Idef.init_image( filename=ref_image_filename)
DIdef.init_image(filename=ref_image_filename)

Martin Genet
committed
DDIdef.init_image(filename=ref_image_filename)
frame_basename = working_folder+"/"+working_basename+"-frame="+str(k_frame).zfill(images_zfill)
if (os.path.exists(frame_basename+".pdf")):
os.remove(frame_basename+".pdf")
file_dat_frame = open(frame_basename+".dat", "w")
file_dat_frame.write("#k_iter B_norm err_res relax dU_norm U_norm err_disp dI_norm err_im err_im_rel\n")
file_pvd_frame = dolfin.File(frame_basename+"_.pvd")
for vtu_filename in glob.glob(frame_basename+"_*.vtu"):
os.remove(vtu_filename)
file_pvd_frame << (U, 0.)
# linear system: matrix

Martin Genet
committed
if (tangent_type.startswith("Idef_old")):
A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})
image_filename = images_folder+"/"+images_basename+"_"+str(k_frame).zfill(images_zfill)+".vti"
Idef.init_image( filename=image_filename)
DIdef.init_image(filename=image_filename)

Martin Genet
committed
DDIdef.init_image(filename=image_filename)
dI_norm = (dolfin.assemble((Idef-Iref)**2 * dX, form_compiler_parameters={'quadrature_degree':images_quadrature})/mesh_volume)**(1./2)
err_im = dI_norm/Iref_norm
file_dat_frame.write(" ".join([str(val) for val in [-2, None, None, None, None, None, None, dI_norm, err_im, None]])+"\n")
dI_norm = (dolfin.assemble((Idef-Iref)**2 * dX, form_compiler_parameters={'quadrature_degree':images_quadrature})/mesh_volume)**(1./2)
err_im = dI_norm/Iref_norm
file_dat_frame.write(" ".join([str(val) for val in [-1, None, None, None, None, None, None, dI_norm, err_im, None]])+"\n")

Martin Genet
committed
if (tangent_type.startswith("Idef")):
A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})
if (relax_type == "aitken"):
if (k_iter == 1):
B_old = B.copy()
elif (k_iter > 1):
B_old[:] = B[:]
B = dolfin.assemble(b, tensor=B, form_compiler_parameters={'quadrature_degree':images_quadrature})
if (relax_type == "constant"):
if (k_iter == 0):
relax = relax_init
elif (relax_type == "aitken"):
if (k_iter == 0):
relax = relax_init
else:
if (k_iter == 1):
dB = B - B_old
elif (k_iter > 1):
dB[:] = B[:] - B_old[:]
relax *= (-1.) * B_old.inner(dB) / dB.inner(dB)
for relax_k in xrange(relax_n_iter):

Martin Genet
committed
U.vector()[:] += (1./relax_n_iter) * dU.vector()[:]
B = dolfin.assemble(b, tensor=B, form_compiler_parameters={'quadrature_degree':images_quadrature})
B_norm_relax[relax_k] = numpy.linalg.norm(B)
print_sci(tab,"B_norm_relax",B_norm_relax[relax_k])

Martin Genet
committed
U.vector()[:] -= dU.vector()[:]
iter_basename = frame_basename+"-iter="+str(k_iter).zfill(3)
open(iter_basename+".dat", "w").write("\n".join([" ".join([str(val) for val in [float(relax_k+1)/relax_n_iter, B_norm_relax[relax_k]]]) for relax_k in xrange(relax_n_iter)]))
os.system("gnuplot -e \"set terminal pdf; set output '"+iter_basename+".pdf'; plot '"+iter_basename+".dat' u 1:2 w l notitle\"")
relax = float(relax_k+1)/relax_n_iter
else:
assert (0), "relax_type must be \"constant\", \"aitken\" or \"manual\". Aborting."

Martin Genet
committed
DU_old.vector()[:] += relax * dU.vector()[:]
U.vector()[:] += relax * dU.vector()[:]
assert (U_norm > 0.), "U_norm = "+str(U_norm)+" <= 0. Aborting."
dI_norm_old = dI_norm
dI_norm = (dolfin.assemble((Idef-Iref)**2 * dX, form_compiler_parameters={'quadrature_degree':images_quadrature})/mesh_volume)**(1./2)
file_dat_frame.write(" ".join([str(val) for val in [k_iter, B_norm, err_res, relax, dU_norm, U_norm, err_disp, dI_norm, err_im, err_im_rel]])+"\n")
if (tol_disp is not None) and (tol_res is not None):
if ((err_disp < tol_disp) or ((dU_norm < tol_disp) and (U_norm < tol_disp))) and (err_res < tol_res):
success = True
elif (tol_disp is not None):
if (err_disp < tol_disp) or ((dU_norm < tol_disp) and (U_norm < tol_disp)):
success = True
elif (tol_res is not None):
if (err_res < tol_res):
success = True
if (success):
if (k_iter == n_iter_max-1):
global_success = False
print_str(tab,"Warning! Nonlinear solver failed to converge… (k_frame = "+str(k_frame)+")")
os.remove(frame_basename+"_.pvd")
file_dat_frame.close()
os.system("gnuplot -e \"set terminal pdf; set output '"+frame_basename+".pdf'; set key box textcolor variable; set grid; set logscale y; set yrange [1e-4:1e1]; plot '"+frame_basename+".dat' u 1:3 pt 1 lw 3 title 'err_res', '' u 1:4 pt 1 lw 3 title 'relax', '' u 1:7 pt 1 lw 3 title 'err_disp', '' using 1:9 pt 1 lw 3 title 'err_im', '' using 1:10 pt 1 lw 3 title 'err_im_rel', "+str(tol_res)+" lt -1 notitle\"")
if not (success) and not (continue_after_fail):
break

Martin Genet
committed
# solution update
U_old.vector()[:] += DU_old.vector()[:]
DU_old.vector()[:] = 0.
if not (success) and not (continue_after_fail):
break