Mentions légales du service

Skip to content
Snippets Groups Projects
fedic.py 17.2 KiB
Newer Older
Martin Genet's avatar
Martin Genet committed
#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

########################################################################

Martin Genet's avatar
Martin Genet committed
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")

Martin Genet's avatar
Martin Genet committed
def fedic(
        working_folder,
        working_basename,
        images_folder,
        images_basename,
        mesh_folder,
        mesh_basename,
Martin Genet's avatar
Martin Genet committed
        images_k_ref=0,
Martin Genet's avatar
Martin Genet committed
        images_zfill=2,
Martin Genet's avatar
Martin Genet committed
        mesh_degree=1,
        penalty=0.9,
        tangent_type="Idef", # Idef, Idef-wHess, Idef_old, Idef_old-wHess, Iref
Martin Genet's avatar
Martin Genet committed
        relax_type="const", # const, aitken, manual
Martin Genet's avatar
Martin Genet committed
        relax_init=1.0,
Martin Genet's avatar
Martin Genet committed
        relax_n_iter=10,
Martin Genet's avatar
Martin Genet committed
        tol_res=None,
        tol_disp=None,
Martin Genet's avatar
Martin Genet committed
        n_iter_max=100,
        print_iterations=0,
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    tab = 0
    print_str(tab,"Loading mesh…")
Martin Genet's avatar
Martin Genet committed
    mesh_filename = mesh_folder+"/"+mesh_basename+".xml"
Martin Genet's avatar
Martin Genet committed
    assert os.path.exists(mesh_filename), "No mesh in "+mesh_filename+". Aborting."
Martin Genet's avatar
Martin Genet committed
    mesh = dolfin.Mesh(mesh_filename)
    dX = dolfin.dx(mesh)
Martin Genet's avatar
Martin Genet committed
    mesh_volume = dolfin.assemble(dolfin.Constant(1)*dX)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Computing quadrature degree for images…")
Martin Genet's avatar
Martin Genet committed
    ref_image_filename = images_folder+"/"+images_basename+"_"+str(images_k_ref).zfill(images_zfill)+".vti"
    if (images_quadrature is None):
        images_quadrature = myFEniCS.compute_quadrature_degree(
            image_filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
            mesh=mesh,
            verbose=0)
    print_var(tab+1,"images_quadrature",images_quadrature)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Loading reference image…")
Martin Genet's avatar
Martin Genet committed
    images_dimension = myVTK.computeImageDimensionality(
Martin Genet's avatar
Martin Genet committed
        verbose=0)
    assert (images_dimension in (2,3)), "images_dimension must be 2 or 3. Aborting."
Martin Genet's avatar
Martin Genet committed
    fe = dolfin.FiniteElement(
        family="Quadrature",
Martin Genet's avatar
Martin Genet committed
        degree=images_quadrature,
        quad_scheme="default")
    ve = dolfin.VectorElement(
        family="Quadrature",
        cell=mesh.ufl_cell(),
Martin Genet's avatar
Martin Genet committed
        degree=images_quadrature,
    te = dolfin.TensorElement(
        family="Quadrature",
        cell=mesh.ufl_cell(),
        degree=images_quadrature,
        quad_scheme="default")
Martin Genet's avatar
Martin Genet committed
    te._quad_scheme = "default"
    for k in xrange(images_dimension**2):
        te.sub_elements()[k]._quad_scheme = "default"
Martin Genet's avatar
Martin Genet committed
    if (images_dimension == 2):
Martin Genet's avatar
Martin Genet committed
        Iref = myFEniCS.ExprIm2(
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
            element=fe)
Martin Genet's avatar
Martin Genet committed
        DIref = myFEniCS.ExprGradIm2(
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
    elif (images_dimension == 3):
Martin Genet's avatar
Martin Genet committed
        Iref = myFEniCS.ExprIm3(
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
            element=fe)
Martin Genet's avatar
Martin Genet committed
        DIref = myFEniCS.ExprGradIm3(
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
    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."
Martin Genet's avatar
Martin Genet committed
    print_var(tab+1,"Iref_norm",Iref_norm)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Defining functions…")
Martin Genet's avatar
Martin Genet committed
    vfs = dolfin.VectorFunctionSpace(
        mesh=mesh,
        family="Lagrange",
        degree=mesh_degree)
Martin Genet's avatar
Martin Genet committed
    U = dolfin.Function(
Martin Genet's avatar
Martin Genet committed
        vfs,
Martin Genet's avatar
Martin Genet committed
        name="displacement")
Martin Genet's avatar
Martin Genet committed
        vfs,
        name="previous displacement")
    DU_old = dolfin.Function(
        vfs,
        name="previous displacement increment")
Martin Genet's avatar
Martin Genet committed
    dU = dolfin.Function(
Martin Genet's avatar
Martin Genet committed
        vfs,
        name="displacement correction")
Martin Genet's avatar
Martin Genet committed
    dU_ = dolfin.TrialFunction(vfs)
    dV_ = dolfin.TestFunction(vfs)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Defining variational forms…")
Martin Genet's avatar
Martin Genet committed
    penalty = dolfin.Constant(penalty)
Martin Genet's avatar
Martin Genet committed
    if (images_dimension == 2):
Martin Genet's avatar
Martin Genet committed
        Idef = myFEniCS.ExprDefIm2(
Martin Genet's avatar
Martin Genet committed
            U=U,
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
            element=fe)
Martin Genet's avatar
Martin Genet committed
        DIdef = myFEniCS.ExprGradDefIm2(
Martin Genet's avatar
Martin Genet committed
            U=U,
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
        DDIdef = myFEniCS.ExprHessDefIm2(
            U=U,
            filename=ref_image_filename,
            element=te)
Martin Genet's avatar
Martin Genet committed
    elif (images_dimension == 3):
Martin Genet's avatar
Martin Genet committed
        Idef = myFEniCS.ExprDefIm3(
Martin Genet's avatar
Martin Genet committed
            U=U,
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
Martin Genet's avatar
Martin Genet committed
            element=fe)
Martin Genet's avatar
Martin Genet committed
        DIdef = myFEniCS.ExprGradDefIm3(
Martin Genet's avatar
Martin Genet committed
            U=U,
Martin Genet's avatar
Martin Genet committed
            filename=ref_image_filename,
        DDIdef = myFEniCS.ExprHessDefIm3(
            U=U,
            filename=ref_image_filename,
            element=te)
Martin Genet's avatar
Martin Genet committed
    if (tangent_type == "Iref"):
        a =     penalty  * dolfin.inner(dolfin.dot(DIref, dU_),
                                        dolfin.dot(DIref, dV_)) * dX
    elif (tangent_type.startswith("Idef")):
Martin Genet's avatar
Martin Genet committed
        a =     penalty  * dolfin.inner(dolfin.dot(DIdef, dU_),
                                        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's avatar
Martin Genet committed
    a +=    (1.-penalty) * dolfin.inner(dolfin.grad(dU_),
                                        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
Martin Genet's avatar
Martin Genet committed
    B0 = dolfin.assemble(b0, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
    B0_norm = numpy.linalg.norm(B0)
Martin Genet's avatar
Martin Genet committed
    assert (B0_norm > 0.), "B0_norm = "+str(B0_norm)+" <= 0. Aborting."
Martin Genet's avatar
Martin Genet committed
    print_var(tab+1,"B0_norm",B0_norm)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    # linear system
    A = None
Martin Genet's avatar
Martin Genet committed
    if (tangent_type == "Iref"):
Martin Genet's avatar
Martin Genet committed
        A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
    B = None
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Checking number of frames…")
    n_frames = len(glob.glob(images_folder+"/"+images_basename+"_"+"?"*images_zfill+".vti"))
Martin Genet's avatar
Martin Genet committed
    #n_frames = 2
    assert (n_frames > 1), "n_frames = "+str(n_frames)+" <= 1. Aborting."
Martin Genet's avatar
Martin Genet committed
    print_var(tab+1,"n_frames",n_frames)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    assert (abs(images_k_ref) < n_frames), "abs(images_k_ref) = "+str(images_k_ref)+" >= n_frames. Aborting."
Martin Genet's avatar
Martin Genet committed
    images_k_ref = images_k_ref%n_frames
Martin Genet's avatar
Martin Genet committed
    print_var(tab+1,"images_k_ref",images_k_ref)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Printing initial solution…")
Martin Genet's avatar
Martin Genet committed
    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))

Martin Genet's avatar
Martin Genet committed
    print_str(tab,"Looping over frames…")
Martin Genet's avatar
Martin Genet committed
    n_iter_tot = 0
Martin Genet's avatar
Martin Genet committed
    global_success = True
Martin Genet's avatar
Martin Genet committed
    for forward_or_backward in ["forward", "backward"]: # not sure if this still works…
Martin Genet's avatar
Martin Genet committed
        print_var(tab,"forward_or_backward",forward_or_backward)
Martin Genet's avatar
Martin Genet committed

        if (forward_or_backward == "forward"):
Martin Genet's avatar
Martin Genet committed
            k_frames = range(images_k_ref+1, n_frames, +1)
Martin Genet's avatar
Martin Genet committed
        elif (forward_or_backward == "backward"):
Martin Genet's avatar
Martin Genet committed
            k_frames = range(images_k_ref-1, -1, -1)
Martin Genet's avatar
Martin Genet committed
        print_var(tab,"k_frames",k_frames)
Martin Genet's avatar
Martin Genet committed

        if (forward_or_backward == "backward"):
            U.vector()[:] = 0.
            U_old.vector()[:] = 0.
            DU_old.vector()[:] = 0.
Martin Genet's avatar
Martin Genet committed
            Idef.init_image( filename=ref_image_filename)
            DIdef.init_image(filename=ref_image_filename)
            DDIdef.init_image(filename=ref_image_filename)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
        tab += 1
Martin Genet's avatar
Martin Genet committed
        for k_frame in k_frames:
Martin Genet's avatar
Martin Genet committed
            print_var(tab-1,"k_frame",k_frame)
Martin Genet's avatar
Martin Genet committed

            if (print_iterations):
Martin Genet's avatar
Martin Genet committed
                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")
Martin Genet's avatar
Martin Genet committed
                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")
Martin Genet's avatar
Martin Genet committed

                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
            if (tangent_type.startswith("Idef_old")):
Martin Genet's avatar
Martin Genet committed
                A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
                #print_var(tab,"A",A.array())
Martin Genet's avatar
Martin Genet committed
                #A_norm = numpy.linalg.norm(A.array())
Martin Genet's avatar
Martin Genet committed
                #print_var(tab,"A_norm",A_norm)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
            print_str(tab,"Loading image and image gradient…")
Martin Genet's avatar
Martin Genet committed
            image_filename = images_folder+"/"+images_basename+"_"+str(k_frame).zfill(images_zfill)+".vti"
Martin Genet's avatar
Martin Genet committed
            Idef.init_image( filename=image_filename)
            DIdef.init_image(filename=image_filename)
            DDIdef.init_image(filename=image_filename)
Martin Genet's avatar
Martin Genet committed

            if (print_iterations):
                U.vector()[:] = 0.
Martin Genet's avatar
Martin Genet committed
                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")
Martin Genet's avatar
Martin Genet committed
                U.vector()[:] = U_old.vector()[:]
Martin Genet's avatar
Martin Genet committed
                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's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
            print_str(tab,"Running registration…")
            tab += 1
Martin Genet's avatar
Martin Genet committed
            k_iter = 0
Martin Genet's avatar
Martin Genet committed
            success = False
Martin Genet's avatar
Martin Genet committed
            while (True):
Martin Genet's avatar
Martin Genet committed
                print_var(tab-1,"k_iter",k_iter)
Martin Genet's avatar
Martin Genet committed
                n_iter_tot += 1

Martin Genet's avatar
Martin Genet committed
                # linear system: matrix
Martin Genet's avatar
Martin Genet committed
                    A = dolfin.assemble(a, tensor=A, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
                    #print_var(tab,"A",A.array())
Martin Genet's avatar
Martin Genet committed
                    #A_norm = numpy.linalg.norm(A.array())
Martin Genet's avatar
Martin Genet committed
                    #print_sci(tab,"A_norm",A_norm)
Martin Genet's avatar
Martin Genet committed

                # linear system: residual
Martin Genet's avatar
Martin Genet committed
                if (relax_type == "aitken"):
                    if (k_iter == 1):
                        B_old = B.copy()
                    elif (k_iter > 1):
                        B_old[:] = B[:]
Martin Genet's avatar
Martin Genet committed
                B = dolfin.assemble(b, tensor=B, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
                #print_var(tab,"B",B.array())
Martin Genet's avatar
Martin Genet committed
                B_norm = numpy.linalg.norm(B.array())
Martin Genet's avatar
Martin Genet committed
                #print_sci(tab,"B_norm",B_norm)
Martin Genet's avatar
Martin Genet committed
                err_res = B_norm/B0_norm
Martin Genet's avatar
Martin Genet committed
                print_sci(tab,"err_res",err_res)
Martin Genet's avatar
Martin Genet committed

                # linear system: solve
Martin Genet's avatar
Martin Genet committed
                dolfin.solve(A, dU.vector(), B)
Martin Genet's avatar
Martin Genet committed
                #print_var(tab,"dU",dU.vector().array())
Martin Genet's avatar
Martin Genet committed

                # relaxation
Martin Genet's avatar
Martin Genet committed
                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)
Martin Genet's avatar
Martin Genet committed
                    print_var(tab,"relax",relax)
Martin Genet's avatar
Martin Genet committed
                elif (relax_type == "manual"):
Martin Genet's avatar
Martin Genet committed
                    B_norm_relax = numpy.empty(relax_n_iter)
                    tab += 1
                    for relax_k in xrange(relax_n_iter):
Martin Genet's avatar
Martin Genet committed
                        print_var(tab-1,"relax_k",relax_k)
                        U.vector()[:] += (1./relax_n_iter) * dU.vector()[:]
Martin Genet's avatar
Martin Genet committed
                        B = dolfin.assemble(b, tensor=B, form_compiler_parameters={'quadrature_degree':images_quadrature})
Martin Genet's avatar
Martin Genet committed
                        B_norm_relax[relax_k] = numpy.linalg.norm(B)
                        print_sci(tab,"B_norm_relax",B_norm_relax[relax_k])
Martin Genet's avatar
Martin Genet committed
                    #print_var(tab,"B_norm_relax",B_norm_relax)
                    if (print_iterations):
Martin Genet's avatar
Martin Genet committed
                        iter_basename = frame_basename+"-iter="+str(k_iter).zfill(3)
Martin Genet's avatar
Martin Genet committed
                        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)]))
Martin Genet's avatar
Martin Genet committed
                        os.system("gnuplot -e \"set terminal pdf; set output '"+iter_basename+".pdf'; plot '"+iter_basename+".dat' u 1:2 w l notitle\"")
Martin Genet's avatar
Martin Genet committed
                    relax_k = numpy.argmin(B_norm_relax)
                    relax = float(relax_k+1)/relax_n_iter
Martin Genet's avatar
Martin Genet committed
                    print_var(tab,"relax",relax)
Martin Genet's avatar
Martin Genet committed
                else:
                    assert (0), "relax_type must be \"constant\", \"aitken\" or \"manual\". Aborting."
Martin Genet's avatar
Martin Genet committed

                # solution update
                DU_old.vector()[:] += relax * dU.vector()[:]
                U.vector()[:] += relax * dU.vector()[:]
Martin Genet's avatar
Martin Genet committed

                if (print_iterations):
Martin Genet's avatar
Martin Genet committed
                    #print_var(tab,"U",U.vector().array())
Martin Genet's avatar
Martin Genet committed
                    file_pvd_frame << (U, float(k_iter+1))
Martin Genet's avatar
Martin Genet committed

                # displacement error
Martin Genet's avatar
Martin Genet committed
                dU_norm = numpy.linalg.norm(dU.vector().array())
Martin Genet's avatar
Martin Genet committed
                print_sci(tab,"dU_norm",dU_norm)
Martin Genet's avatar
Martin Genet committed
                U_norm = numpy.linalg.norm(U.vector().array())
Martin Genet's avatar
Martin Genet committed
                assert (U_norm > 0.), "U_norm = "+str(U_norm)+" <= 0. Aborting."
Martin Genet's avatar
Martin Genet committed
                print_sci(tab,"U_norm",U_norm)
Martin Genet's avatar
Martin Genet committed
                err_disp = dU_norm/U_norm
Martin Genet's avatar
Martin Genet committed
                print_sci(tab,"err_disp",err_disp)
Martin Genet's avatar
Martin Genet committed

                # image error
Martin Genet's avatar
Martin Genet committed
                if (k_iter > 0):
Martin Genet's avatar
Martin Genet committed
                    dI_norm_old = dI_norm
                dI_norm = (dolfin.assemble((Idef-Iref)**2 * dX, form_compiler_parameters={'quadrature_degree':images_quadrature})/mesh_volume)**(1./2)
Martin Genet's avatar
Martin Genet committed
                #print_sci(tab,"dI_norm",dI_norm)
Martin Genet's avatar
Martin Genet committed
                err_im = dI_norm/Iref_norm
Martin Genet's avatar
Martin Genet committed
                print_sci(tab,"err_im",err_im)
Martin Genet's avatar
Martin Genet committed
                if (k_iter == 0):
                    err_im_rel = 1.
                else:
Martin Genet's avatar
Martin Genet committed
                    err_im_rel = abs(dI_norm-dI_norm_old)/dI_norm_old
Martin Genet's avatar
Martin Genet committed
                    print_sci(tab,"err_im_rel",err_im_rel)
Martin Genet's avatar
Martin Genet committed

                if (print_iterations):
Martin Genet's avatar
Martin Genet committed
                    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")
Martin Genet's avatar
Martin Genet committed

                # exit test
Martin Genet's avatar
Martin Genet committed
                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):
Martin Genet's avatar
Martin Genet committed
                    print_str(tab,"Nonlinear solver converged…")
Martin Genet's avatar
Martin Genet committed
                    break

Martin Genet's avatar
Martin Genet committed
                if (k_iter == n_iter_max-1):
                    global_success = False
Martin Genet's avatar
Martin Genet committed
                    print_str(tab,"Warning! Nonlinear solver failed to converge… (k_frame = "+str(k_frame)+")")
Martin Genet's avatar
Martin Genet committed
                    break

                # increment counter
                k_iter += 1

Martin Genet's avatar
Martin Genet committed
            tab -= 1

Martin Genet's avatar
Martin Genet committed
            if (print_iterations):
Martin Genet's avatar
Martin Genet committed
                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\"")
Martin Genet's avatar
Martin Genet committed

            if not (success) and not (continue_after_fail):
                break

            # solution update
            U_old.vector()[:] += DU_old.vector()[:]
            DU_old.vector()[:] = 0.

Martin Genet's avatar
Martin Genet committed
            print_str(tab,"Printing solution…")
Martin Genet's avatar
Martin Genet committed
            file_pvd << (U, float(k_frame))
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
        tab -= 1

Martin Genet's avatar
Martin Genet committed
        if not (success) and not (continue_after_fail):
            break

Martin Genet's avatar
Martin Genet committed
    print_var(tab,"n_iter_tot",n_iter_tot)
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    os.remove(pvd_basename+".pvd")
Martin Genet's avatar
Martin Genet committed

Martin Genet's avatar
Martin Genet committed
    return global_success