This page was generated from unit-3.8-nonlmin/nonlmin.ipynb.
3.8 Nonlinear minimization problems¶
We consider problems of the form
[1]:
from ngsolve import *
import netgen.gui
Scalar minimization problems¶
As a first example we take \(V = H^1_0\) and
The minimization is equivalent to solving the nonlinear PDE:
We solve the PDE with a Newton iteration.
[2]:
from netgen.geom2d import unit_square
mesh = Mesh (unit_square.GenerateMesh(maxh=0.2))
V = H1(mesh, order=4, dirichlet=[1,2,3,4])
u = V.TrialFunction()
Generate Mesh from spline geometry
Boundary mesh done, np = 20
CalcLocalH: 20 Points 0 Elements 0 Surface Elements
Meshing domain 1 / 1
load internal triangle rules
Surface meshing done
Edgeswapping, topological
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Update mesh topology
Update clusters
To solve the problem we use the Variation
integrator. Based on the symbolic description of the energy functional, it is able to
evaluate the energy functional (
Energy
)
compute the Gateau derivative for a given \(u\) (
Apply
):
compute the second derivative (
AssembleLinearization
)
[3]:
a = BilinearForm (V, symmetric=True)
a += Variation ( (grad(u)*grad(u) + u**4-u) * dx)
Equivalent to:
a += (2 * grad(u) * grad(v) + 4*u*u*u*v - 1 * v)*dx
(which has the same form as the problems in the nonlinear example)
We recall the Newton iteration (cf. unit-3.7 ) we make the loop:
Given an initial guess \(u^0\)
loop over \(i=0,..\) until convergence:
Compute linearization: $A u^i + \delta `A(u^i) :nbsphinx-math:Delta `u^{i} = 0 $:
\(f^i = A u^i\)
\(B^i = \delta A(u^i)\)
Solve \(B^i \Delta u^i = -f^i\)
Update \(u^{i+1} = u^i + \Delta u^{i}\)
Evaluate stopping criteria
Evaluate \(E(u^{i+1})\)
As a stopping criteria we take \(\langle A u^i,\Delta u^i \rangle = \langle A u^i, A u^i \rangle_{(B^i)^{-1}}< \varepsilon\).
[4]:
def SolveNonlinearMinProblem(a,gfu,tol=1e-13,maxits=25):
res = gfu.vec.CreateVector()
du = gfu.vec.CreateVector()
for it in range(maxits):
print ("Newton iteration {:3}".format(it),end="")
print ("energy = {:16}".format(a.Energy(gfu.vec)),end="")
#solve linearized problem:
a.Apply (gfu.vec, res)
a.AssembleLinearization (gfu.vec)
inv = a.mat.Inverse(V.FreeDofs())
du.data = inv * res
#update iteration
gfu.vec.data -= du
#stopping criteria
stopcritval = sqrt(abs(InnerProduct(du,res)))
print ("<A u",it,", A u",it,">_{-1}^0.5 = ", stopcritval)
if stopcritval < tol:
break
Redraw(blocking=True)
[5]:
gfu = GridFunction (V)
gfu.vec[:] = 0
Draw(gfu,mesh,"u")
SolveNonlinearMinProblem(a,gfu)
print ("energy = ", a.Energy(gfu.vec))
Newton iteration 0energy = 0.0Assemble linearization
assemble VOL element 54/54
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.13255957337312613
Newton iteration 1energy = -0.008785676900152117Assemble linearization
assemble VOL element 54/54
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.1107598189167146e-05
Newton iteration 2energy = -0.008785676961841486Assemble linearization
assemble VOL element 54/54
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.807517603399629e-13
Newton iteration 3energy = -0.008785676961841488Assemble linearization
assemble VOL element 54/54
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 3.3304484794486327e-17
energy = -0.008785676961841488
Again, a Newton for minimization is shipped with NGSolve:
[6]:
from ngsolve.solvers import *
gfu.vec[:] = 0
NewtonMinimization(a,gfu)
Redraw()
Newton iteration 0
Energy: 0.0
Assemble linearization
assemble VOL element 54/54
call umfpack ... done
err = 0.13255957337312615
Newton iteration 1
Energy: -0.008785676900152117
Assemble linearization
assemble VOL element 54/54
call umfpack ... done
err = 1.1107598189171795e-05
Newton iteration 2
Energy: -0.008785676961841486
Assemble linearization
assemble VOL element 54/54
call umfpack ... done
err = 2.80747021639635e-13
Nonlinear elasticity¶
We consider a beam which is fixed on one side and is subject to gravity only. We assume a Neo-Hookean hyperelastic material. The model is a nonlinear minimization problem with
where \(\mu\) and \(\lambda\) are the Lamé parameters and \(F = I + D v\) where \(v: \Omega \to \mathbb{R}^2\) is the sought for displacement.
[7]:
import netgen.geom2d as geom2d
from ngsolve import *
geo = geom2d.SplineGeometry()
pnums = [ geo.AddPoint (x,y,maxh=0.01) for x,y in [(0,0), (1,0), (1,0.1), (0,0.1)] ]
for p1,p2,bc in [(0,1,"bot"), (1,2,"right"), (2,3,"top"), (3,0,"left")]:
geo.Append(["line", pnums[p1], pnums[p2]], bc=bc)
mesh = Mesh(geo.GenerateMesh(maxh=0.05))
Generate Mesh from spline geometry
Boundary mesh done, np = 68
CalcLocalH: 68 Points 0 Elements 0 Surface Elements
Meshing domain 1 / 1
Surface meshing done
Edgeswapping, topological
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Update mesh topology
Update clusters
[8]:
# E module and poisson number:
E, nu = 210, 0.2
# Lamé constants:
mu = E / 2 / (1+nu)
lam = E * nu / ((1+nu)*(1-2*nu))
V = VectorH1(mesh, order=2, dirichlet="left")
u = V.TrialFunction()
#gravity:
force = CoefficientFunction( (0,-1) )
[9]:
def Pow(a, b):
return exp (log(a)*b)
def NeoHook (C):
return 0.5 * mu * (Trace(C-I) + 2*mu/lam * Pow(Det(C), -lam/2/mu) - 1)
I = Id(mesh.dim)
F = I + Grad(u)
C = F.trans * F
factor = Parameter(1.0)
a = BilinearForm(V, symmetric=True)
a += Variation( NeoHook (C).Compile() * dx
-factor * (InnerProduct(force,u) ).Compile() * dx)
Compiled CF:
N5ngfem27IdentityCoefficientFunctionE
N5ngfem13ProxyFunctionE
N5ngfem13cl_BinaryOpCFINS_11GenericPlusEEE
N5ngfem28TransposeCoefficientFunctionE
N5ngfem29MultMatMatCoefficientFunctionE
N5ngfem13cl_BinaryOpCFINS_12GenericMinusEEE
N5ngfem24TraceCoefficientFunctionE
N5ngfem30DeterminantCoefficientFunctionILi2EEE
N5ngfem12cl_UnaryOpCFI10GenericLogEE
N5ngfem24ScaleCoefficientFunctionE
N5ngfem12cl_UnaryOpCFI10GenericExpEE
N5ngfem24ScaleCoefficientFunctionE
N5ngfem13cl_BinaryOpCFINS_11GenericPlusEEE
N5ngfem27ConstantCoefficientFunctionE
N5ngfem13cl_BinaryOpCFINS_12GenericMinusEEE
N5ngfem24ScaleCoefficientFunctionE
inputs =
0:
1:
2: 0 1
3: 2
4: 3 2
5: 4 0
6: 5
7: 4
8: 7
9: 8
10: 9
11: 10
12: 6 11
13:
14: 12 13
15: 14
Compiled CF:
N5ngfem27ConstantCoefficientFunctionE
N5ngfem27ConstantCoefficientFunctionE
N5ngfem28VectorialCoefficientFunctionE
N5ngfem12cl_UnaryOpCFI15GenericIdentityEE
N5ngfem13ProxyFunctionE
N5ngfem31T_MultVecVecCoefficientFunctionILi2EEE
inputs =
0:
1:
2: 0 1
3: 2
4:
5: 3 4
We want to solve the minimization problem for \(\gamma = 5\). Due to the high nonlinearity in the problem, the Newton iteration will not convergence with any initial guess. We approach the case \(\gamma = 5\) by solving problems with \(\gamma = i/10\) for \(i=1,..,50\) and taking the solution of the previous problem as an initial guess.
[10]:
gfu = GridFunction(V)
gfu.vec[:] = 0
Draw (gfu, mesh, "u")
SetVisualization (deformation=True)
res = gfu.vec.CreateVector()
du = gfu.vec.CreateVector()
for loadstep in range(50):
print ("loadstep", loadstep)
factor.Set ((loadstep+1)/10)
SolveNonlinearMinProblem(a,gfu)
Redraw()
loadstep 0
Newton iteration 0energy = 8.749999999999998Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.01667849444523386
Newton iteration 1energy = 8.750132591453418Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.023333560516275564
Newton iteration 2energy = 8.74986114895739Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00010232534273250714
Newton iteration 3energy = 8.749861143722152Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.1484913353753e-08
Newton iteration 4energy = 8.749861143722145Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.734694947345281e-13
Newton iteration 5energy = 8.749861143722145Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 9.247273008118622e-16
loadstep 1
Newton iteration 0energy = 8.74958388866222Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.016596172790229985
Newton iteration 1energy = 8.749710096254224Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.022958523347670023
Newton iteration 2energy = 8.749447295506688Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00014215759220200522
Newton iteration 3energy = 8.749447285508145Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.325574122106694e-06
Newton iteration 4energy = 8.749447285505438Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 6.329020240483118e-10
Newton iteration 5energy = 8.749447285505443Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 1.6967201531172987e-15
loadstep 2
Newton iteration 0energy = 8.748898140386652Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.016356983393359056
Newton iteration 1energy = 8.749005259290692Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.021892894044873473
Newton iteration 2energy = 8.748766252622183Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.000208222077261969
Newton iteration 3energy = 8.748766231266895Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.890915769550513e-06
Newton iteration 4energy = 8.748766231254931Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.665249442467225e-09
Newton iteration 5energy = 8.748766231254935Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 2.2360660986729218e-15
loadstep 3
Newton iteration 0energy = 8.747955292062146Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0159820733778304
Newton iteration 1energy = 8.74803542180493Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.020292399713242315
Newton iteration 2energy = 8.747830036432735Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00026380399717611486
Newton iteration 3energy = 8.747830002128063Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 7.24731483392178e-06
Newton iteration 4energy = 8.747830002101788Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 4.744152101219956e-09
Newton iteration 5energy = 8.747830002101791Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 3.834758384860599e-15
loadstep 4
Newton iteration 0energy = 8.746771010509391Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.015501067150346784
Newton iteration 1energy = 8.746821856285143Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.01835991696222946
Newton iteration 2energy = 8.746653679174742Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00029906276685253894
Newton iteration 3energy = 8.746653634982097Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 8.720360402018377e-06
Newton iteration 4energy = 8.746653634944062Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 5.243058230169497e-09
Newton iteration 5energy = 8.746653634944062Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 4.57624018420641e-15
loadstep 5
Newton iteration 0energy = 8.745362710869802Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.01494623240146039
Newton iteration 1energy = 8.745386462940143Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.016292996834189524
Newton iteration 2energy = 8.74525397874261Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00031320373730301293
Newton iteration 3energy = 8.745253930142317Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 9.042918405812291e-06
Newton iteration 4energy = 8.745253930101422Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 4.2170747228818445e-09
Newton iteration 5energy = 8.74525393010142Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 4.4934195509996876e-15
loadstep 6
Newton iteration 0energy = 8.743748369013607Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.014347713355602996
Newton iteration 1energy = 8.743749760522784Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.014249104212014871
Newton iteration 2energy = 8.74364839737191Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00030934606309694624
Newton iteration 3energy = 8.743648349846454Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 8.384888859530472e-06
Newton iteration 4energy = 8.743648349811291Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.6940148501867647e-09
Newton iteration 5energy = 8.7436483498113Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 4.577084252820437e-15
loadstep 7
Newton iteration 0energy = 8.741945639262823Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.013730721530226796
Newton iteration 1energy = 8.741930180324744Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.012333575792311013
Newton iteration 2energy = 8.741854213024633Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.0002921321706594058
Newton iteration 3energy = 8.74185417055718Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 7.130071890016913e-06
Newton iteration 4energy = 8.741854170531763Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.4434881033484502e-09
Newton iteration 5energy = 8.74185417053177Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 5.262958683534736e-15
loadstep 8
Newton iteration 0energy = 8.739971280231286Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.013114572255962758
Newton iteration 1energy = 8.73994411367146Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.010604017031859161
Newton iteration 2energy = 8.739887940284762Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00026632144941452456
Newton iteration 3energy = 8.73988790493644Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.662818925120992e-06
Newton iteration 4energy = 8.739887904920405Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 6.720247861727614e-10
Newton iteration 5energy = 8.739887904920398Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 5.5915825831605006e-15
loadstep 9
Newton iteration 0energy = 8.73784083403977Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0125129674198164
Newton iteration 1energy = 8.737806204166365Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.009082089187694535
Newton iteration 2energy = 8.73776498517593Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00023604345317251825
Newton iteration 3energy = 8.73776495737735Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.2585795846192444e-06
Newton iteration 4energy = 8.737764957368277Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.778882710277497e-10
Newton iteration 5energy = 8.73776495736828Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 6.192838256475658e-15
loadstep 10
Newton iteration 0energy = 8.735568488681253Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.011934895047581391
Newton iteration 1energy = 8.735529626298675Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.007766086303887108
Newton iteration 2energy = 8.735499478251771Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00020450110958987553
Newton iteration 3energy = 8.735499457369832Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 3.0649641069813007e-06
Newton iteration 4energy = 8.73549945736514Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0317259419213303e-10
Newton iteration 5energy = 8.735499457365144Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 6.1435177694761e-15
loadstep 11
Newton iteration 0energy = 8.733167060537248Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.011385697188469452
Newton iteration 1energy = 8.733126283518718Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.006641150145740644
Newton iteration 2energy = 8.733104230848486Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00017394439104879262
Newton iteration 3energy = 8.733104215732606Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.129411853086759e-06
Newton iteration 4energy = 8.733104215730345Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 3.432558772877624e-11
Newton iteration 5energy = 8.733104215730338Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 8.094250286114186e-15
loadstep 12
Newton iteration 0energy = 8.730648048952677Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.01086806034480337
Newton iteration 1energy = 8.730606940950912Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.005686383297879251
Newton iteration 2energy = 8.730590769249666Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00014578425109884895
Newton iteration 3energy = 8.730590758628075Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.4384336718370366e-06
Newton iteration 4energy = 8.730590758627047Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.987280236896591e-12
Newton iteration 5energy = 8.730590758627043Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 7.875973513059426e-15
loadstep 13
Newton iteration 0energy = 8.72802172950889Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.01038282956646614
Newton iteration 1energy = 8.727981322866317Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.004879270983755723
Newton iteration 2energy = 8.727969413471362Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 0.00012076411719468836
Newton iteration 3energy = 8.727969406181067Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 9.505522832199866e-07
Newton iteration 4energy = 8.727969406180618Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.327209811059151e-12
Newton iteration 5energy = 8.727969406180616Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 8.366395854739112e-15
loadstep 14
Newton iteration 0energy = 8.725297264364867Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.009929630827329049
Newton iteration 1energy = 8.725258196156323Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.004198149297602768
Newton iteration 2energy = 8.72524937798614Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.914019138175112e-05
Newton iteration 3energy = 8.725249373072183Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 6.177438537851986e-07
Newton iteration 4energy = 8.725249373071996Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.711757828663729e-13
Newton iteration 5energy = 8.72524937307199Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 8.086506222001764e-15
loadstep 15
Newton iteration 0energy = 8.722482816430372Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.009507324233198547
Newton iteration 1energy = 8.72244544987052Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0036234023936600923
Newton iteration 2energy = 8.722438879876965Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.084368163731423e-05
Newton iteration 3energy = 8.722438876609113Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 3.966066049271973e-07
Newton iteration 4energy = 8.722438876609038Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.5278000509753718e-13
Newton iteration 5energy = 8.722438876609038Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 8.08115538661459e-15
loadstep 16
Newton iteration 0energy = 8.719585659767372Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.009114322388199479
Newton iteration 1energy = 8.71955017364448Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00313789936072171
Newton iteration 2energy = 8.719545245704236Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 6.561232996881842e-05
Newton iteration 3energy = 8.71954524355164Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.52533093944595e-07
Newton iteration 4energy = 8.719545243551611Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.5033040334992298e-13
Newton iteration 5energy = 8.719545243551606Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 5 , A u 5 >_{-1}^0.5 = 8.051020566432653e-15
loadstep 17
Newton iteration 0energy = 8.71661228222845Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.008748807494135346
Newton iteration 1energy = 8.716578734485331Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0027270119708742158
Newton iteration 2energy = 8.716575012244606Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.308808014174092e-05
Newton iteration 3energy = 8.716575010835323Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.599966209600588e-07
Newton iteration 4energy = 8.716575010835307Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.245681199048317e-14
loadstep 18
Newton iteration 0energy = 8.713568478570464Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.008408875280311978
Newton iteration 1energy = 8.713536850394119Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0023784245288753533
Newton iteration 2energy = 8.713534018737166Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.288364628582515e-05
Newton iteration 3energy = 8.71353401781759Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.011410754526628e-07
Newton iteration 4energy = 8.71353401781757Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 4.8832450515969854e-14
loadstep 19
Newton iteration 0energy = 8.710459433603457Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.008092627477923432
Newton iteration 1energy = 8.710429659421148Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.002081858695703015
Newton iteration 2energy = 8.710427489790115Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.462374278328321e-05
Newton iteration 3energy = 8.71042748919065Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 6.393601985494953e-08
Newton iteration 4energy = 8.710427489190652Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.505219063009473e-14
loadstep 20
Newton iteration 0energy = 8.707289795678951Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.007798228855193676
Newton iteration 1energy = 8.707261783232546Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.001828780371476849
Newton iteration 2energy = 8.707260108985894Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.796754346607161e-05
Newton iteration 3energy = 8.707260108594774Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.0490580515779663e-08
Newton iteration 4energy = 8.707260108594774Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.4094429423583045e-14
loadstep 21
Newton iteration 0energy = 8.704063741206715Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.007523940275902371
Newton iteration 1energy = 8.704037384768432Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0016121223324776503
Newton iteration 2energy = 8.704036083705963Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.2618373003579114e-05
Newton iteration 3energy = 8.70403608345015Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.57264200426525e-08
Newton iteration 4energy = 8.704036083450148Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.025503816822179e-14
loadstep 22
Newton iteration 0energy = 8.700785031058096Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.007268135817663741
Newton iteration 1energy = 8.700760219966702Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0014260372524824256
Newton iteration 2energy = 8.700759201928879Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.8325434222126603e-05
Newton iteration 3energy = 8.700759201760949Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.6417609086959264e-08
Newton iteration 4energy = 8.70075920176095Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.90254098112901e-15
loadstep 23
Newton iteration 0energy = 8.697457059755887Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.007029309502644591
Newton iteration 1energy = 8.697433683793736Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0012656853969031537
Newton iteration 2energy = 8.69743288184044Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.4881058448227105e-05
Newton iteration 3energy = 8.69743288172972Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.05321833182215e-08
Newton iteration 4energy = 8.697432881729718Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.588291840221449e-15
loadstep 24
Newton iteration 0energy = 8.694082898322856Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.006806075429231539
Newton iteration 1energy = 8.694060850980257Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.001127056033072989
Newton iteration 2energy = 8.69406021509277Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.2115810532600391e-05
Newton iteration 3energy = 8.694060215019382Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 6.796446181591913e-09
Newton iteration 4energy = 8.694060215019384Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0653141855886985e-14
loadstep 25
Newton iteration 0energy = 8.690665331596561Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.006597163862847386
Newton iteration 1energy = 8.690644511937402Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.001006819235980881
Newton iteration 2energy = 8.69064400450148Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.892890798941409e-06
Newton iteration 3energy = 8.690644004452542Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.413672136741346e-09
Newton iteration 4energy = 8.690644004452544Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 8.706654914525168e-15
loadstep 26
Newton iteration 0energy = 8.687206890738581Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.00640141499389159
Newton iteration 1energy = 8.68718720435108Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0009022039607243906
Newton iteration 2energy = 8.687186796900923Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.102650665520232e-06
Newton iteration 3energy = 8.68718679686809Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.8854231418118854e-09
Newton iteration 4energy = 8.687186796868094Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.50793693421982e-15
loadstep 27
Newton iteration 0energy = 8.683709881582512Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.006217771488356258
Newton iteration 1energy = 8.683691240941258Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0008108982398810496
Newton iteration 2energy = 8.683690911799431Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 6.6576303877916165e-06
Newton iteration 3energy = 8.683690911777271Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.8993156161478608e-09
Newton iteration 4energy = 8.683690911777267Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.362268095886705e-15
loadstep 28
Newton iteration 0energy = 8.680176409385007Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.006045270560225572
Newton iteration 1energy = 8.680158733842598Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0007309677194524074
Newton iteration 2energy = 8.680158466399725Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.488278554618818e-06
Newton iteration 3energy = 8.680158466384661Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.2589535377386467e-09
Newton iteration 4energy = 8.680158466384663Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.90942431883319e-15
loadstep 29
Newton iteration 0energy = 8.676608400469208Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005883036026966293
Newton iteration 1energy = 8.676591616023208Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0006607892405581196
Newton iteration 2energy = 8.676591397476551Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.539372421135213e-06
Newton iteration 3energy = 8.676591397466245Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 8.403590576691722e-10
Newton iteration 4energy = 8.676591397466247Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.866077327214238e-15
loadstep 30
Newton iteration 0energy = 8.673007621183567Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005730270629910569
Newton iteration 1energy = 8.67299166011491Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0005989966873394284
Newton iteration 2energy = 8.67299148053809Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.7670865979156124e-06
Newton iteration 3energy = 8.672991480530989Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.648759494872811e-10
Newton iteration 4energy = 8.672991480530989Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0050875893853173e-14
loadstep 31
Newton iteration 0energy = 8.66937569453993Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0055862487812907875
Newton iteration 1energy = 8.669360494985764Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0005444367972390969
Newton iteration 2energy = 8.66936034663853Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.136627189130437e-06
Newton iteration 3energy = 8.669360346633617Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 3.8233940124142637e-10
Newton iteration 4energy = 8.669360346633615Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1734935299802074e-14
loadstep 32
Newton iteration 0energy = 8.665714114843603Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005450309821622411
Newton iteration 1energy = 8.665699620344926Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0004961330477172054
Newton iteration 2energy = 8.665699497158188Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.620341391709445e-06
Newton iteration 3energy = 8.665699497154765Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.6056570946754966e-10
Newton iteration 4energy = 8.665699497154764Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.217201318445093e-14
loadstep 33
Newton iteration 0energy = 8.662024260583653Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005321851819893844
Newton iteration 1energy = 8.662010419633267Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0004532560862173065
Newton iteration 2energy = 8.662010316822602Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.1962173093012577e-06
Newton iteration 3energy = 8.662010316820195Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.787715409258059e-10
Newton iteration 4energy = 8.662010316820183Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.067211648222395e-14
loadstep 34
Newton iteration 0energy = 8.658307405813543Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005200325917091465
Newton iteration 1energy = 8.658294171419934Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0004150994634587661
Newton iteration 2energy = 8.65829408519383Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.8466987982871333e-06
Newton iteration 3energy = 8.65829408519212Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.2346713844118924e-10
Newton iteration 4energy = 8.658294085192125Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 9.32755969194597e-15
loadstep 35
Newton iteration 0energy = 8.654564730219906Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.005085231193936819
Newton iteration 1energy = 8.654552059495636Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00038105966948391123
Newton iteration 2energy = 8.654551986834138Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.557751764245773e-06
Newton iteration 3energy = 8.654551986832924Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 8.582341612129994e-11
Newton iteration 4energy = 8.654551986832924Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.2058191498589807e-14
loadstep 36
Newton iteration 0energy = 8.650797328048808Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004976110032586554
Newton iteration 1energy = 8.65078518182817Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00035061966629611104
Newton iteration 2energy = 8.650785120313985Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.3181296369092865e-06
Newton iteration 3energy = 8.650785120313117Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 6.003349278006359e-11
Newton iteration 4energy = 8.650785120313119Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0325850670695882e-14
loadstep 37
Newton iteration 0energy = 8.647006216035676Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004872543936238474
Newton iteration 1energy = 8.646994558523108Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00032333526738206875
Newton iteration 2energy = 8.646994506212032Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.1187957572450256e-06
Newton iteration 3energy = 8.646994506211405Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.2252692019247613e-11
Newton iteration 4energy = 8.646994506211412Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0729852040909942e-14
loadstep 38
Newton iteration 0energy = 8.64319234046425Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004774149768405117
Newton iteration 1energy = 8.64318113891363Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00029882384028847404
Newton iteration 2energy = 8.643181094234647Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.52469007074102e-07
Newton iteration 3energy = 8.643181094234198Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.99157638633773e-11
Newton iteration 4energy = 8.643181094234194Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.231791869358724e-14
loadstep 39
Newton iteration 0energy = 8.639356583463016Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004680576373940036
Newton iteration 1energy = 8.639345807886512Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0002767549095867613
Newton iteration 2energy = 8.639345769564418Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.132660428051996e-07
Newton iteration 3energy = 8.639345769564084Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.1305856722228868e-11
Newton iteration 4energy = 8.639345769564079Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1818176358974371e-14
loadstep 40
Newton iteration 0energy = 8.635499768632432Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004591501544822213
Newton iteration 1energy = 8.635489391537332Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0002568423180182777
Newton iteration 2energy = 8.635489358532457Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 6.96419178009762e-07
Newton iteration 3energy = 8.635489358532212Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.5262192556902303e-11
Newton iteration 4energy = 8.635489358532213Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1024939938001828e-14
loadstep 41
Newton iteration 0energy = 8.631622666083985Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004506629296485266
Newton iteration 1energy = 8.631612662235378Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0002388376690438732
Newton iteration 2energy = 8.631612633696458Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.980535273693957e-07
Newton iteration 3energy = 8.631612633696278Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.099135182257061e-11
Newton iteration 4energy = 8.631612633696287Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1420765638026738e-14
loadstep 42
Newton iteration 0energy = 8.627725996961052Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0044256874227469566
Newton iteration 1energy = 8.627716343168219Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00022252482582335286
Newton iteration 2energy = 8.627716318395358Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.150106027527706e-07
Newton iteration 3energy = 8.627716318395228Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 7.959435264709924e-12
Newton iteration 4energy = 8.627716318395228Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1578943244195082e-14
loadstep 43
Newton iteration 0energy = 8.62381043750256Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0043484253001886426
Newton iteration 1energy = 8.62380111242687Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00020771528351786156
Newton iteration 2energy = 8.623801090842262Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.447083520833389e-07
Newton iteration 3energy = 8.623801090842166Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.793345816132156e-12
Newton iteration 4energy = 8.623801090842163Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.0945054837429887e-14
loadstep 44
Newton iteration 0energy = 8.619876622702083Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004274611915556383
Newton iteration 1energy = 8.619867606684437Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0001942442655680472
Newton iteration 2energy = 8.619867587809216Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.850298436002284e-07
Newton iteration 3energy = 8.619867587809129Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.240958949701815e-12
Newton iteration 4energy = 8.619867587809136Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1418242765222334e-14
loadstep 45
Newton iteration 0energy = 8.615925149608909Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004204034092581162
Newton iteration 1energy = 8.615916424514648Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00018196742192939653
Newton iteration 2energy = 8.615916407950385Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.3423451071494996e-07
Newton iteration 3energy = 8.61591640795033Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 3.1185096160350177e-12
Newton iteration 4energy = 8.615916407950332Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.110429334265253e-14
loadstep 46
Newton iteration 0energy = 8.611956580310858Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004136494896624677
Newton iteration 1energy = 8.611948129390502Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0001707580288323152
Newton iteration 2energy = 8.611948114804498Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.9088719186477544e-07
Newton iteration 3energy = 8.611948114804452Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.305393316626788e-12
Newton iteration 4energy = 8.611948114804456Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.2157063268037797e-14
loadstep 47
Newton iteration 0energy = 8.607971444634586Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004071812198364899
Newton iteration 1energy = 8.607963252398571Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0001605046079768612
Newton iteration 2energy = 8.607963239511943Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.5380125196407897e-07
Newton iteration 3energy = 8.60796323951191Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.710252507370915e-12
Newton iteration 4energy = 8.607963239511916Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.2290069453583734e-14
loadstep 48
Newton iteration 0energy = 8.603970242593698Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0040098173792995325
Newton iteration 1energy = 8.603962294699814Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00015110889702384532
Newton iteration 2energy = 8.60396228327801Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.219928560793886e-07
Newton iteration 3energy = 8.603962283277987Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 1.2781909165814276e-12
Newton iteration 4energy = 8.603962283277989Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.2774778602252051e-14
loadstep 49
Newton iteration 0energy = 8.599953446612274Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.003950354164065159
Newton iteration 1energy = 8.599945729764308Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00014248411538807357
Newton iteration 2energy = 8.599945719609346Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.9464412525932455e-07
Newton iteration 3energy = 8.599945719609325Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 9.577059269412203e-13
Newton iteration 4energy = 8.599945719609325Assemble linearization
assemble VOL element 170/170
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 1.1386571464267781e-14
Allen-Cahn equation¶
The Allen-Cahn equations describe the process of phase separation and is the (\(L^2\)) gradient-flow equation to the energy
i.e. the solution to the Allen-Cahn equation solves
The quantity \(u\) is an indicator for a phase where \(-1\) refers to one phase and \(1\) to another phase.
The equation has two driving forces:
\(u\) is pulled into one of the two minima (\(-1\) and \(1\)) of the nonlinear term \(u^2(1-u^2)\) (separation of the phases)
the diffusion term scaled with \(\varepsilon\) enforces a smooth transition between the two phases. \(\varepsilon\) determines the size of the transition layer
We use the “SymbolicEnergy” feature to formulate the energy minimization problem and combine it with an implicit Euler discretization:
which we can interpreted as a nonlinear minimization problem again with the energy
To solve the nonlinear equation at every time step we again rely on Newton’s method.
[11]:
from ngsolve import *
from netgen.geom2d import *
periodic = SplineGeometry()
pnts = [ (0,0), (1,0), (1,1), (0,1) ]
pnums = [periodic.AppendPoint(*p) for p in pnts]
lright = periodic.Append ( ["line", pnums[0], pnums[1]],bc="periodic")
btop = periodic.Append ( ["line", pnums[1], pnums[2]], bc="periodic")
periodic.Append ( ["line", pnums[3], pnums[2]], leftdomain=0, rightdomain=1, copy=lright, bc="periodic")
periodic.Append ( ["line", pnums[0], pnums[3]], leftdomain=0, rightdomain=1, copy=btop, bc="periodic")
mesh = Mesh (periodic.GenerateMesh(maxh=0.2))
V = Periodic(H1(mesh, order=4, dirichlet=[]))
u = V.TrialFunction()
eps = 4e-3
dt = 1e-1
gfu = GridFunction(V)
gfuold = GridFunction(V)
a = BilinearForm (V, symmetric=False)
a += Variation( (eps/2*grad(u)*grad(u) + ((1-u*u)*(1-u*u))
+ 0.5/dt*(u-gfuold)*(u-gfuold)) * dx)
Generate Mesh from spline geometry
Copy edge, from 1 to 3
Copy edge, from 2 to 4
Boundary mesh done, np = 20
CalcLocalH: 20 Points 0 Elements 0 Surface Elements
Meshing domain 1 / 1
Surface meshing done
Edgeswapping, topological
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Update mesh topology
Update clusters
[12]:
from math import pi
gfu = GridFunction(V)
gfu.Set(sin(2*pi*x))
#gfu.Set(sin(1e8*x)) #<- essentially a random function
Draw(gfu,mesh,"u")
SetVisualization (deformation=True)
t = 0
setvalues element 58/58
[13]:
for timestep in range(50):
gfuold.vec.data = gfu.vec
SolveNonlinearMinProblem(a,gfu)
Redraw()
t += dt
print("t = ", t)
Newton iteration 0energy = 0.4145019207496362Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.28092704638147636
Newton iteration 1energy = 0.37689785802772957Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.019174370471939796
Newton iteration 2energy = 0.3767134210404821Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.628937134596328e-05
Newton iteration 3energy = 0.3767134164045804Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.505456303796458e-09
Newton iteration 4energy = 0.3767134164045804Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 2.8167247591431503e-16
t = 0.1
Newton iteration 0energy = 0.34388128652385935Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.2360393474078292
Newton iteration 1energy = 0.31710944329874896Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.013442999343193951
Newton iteration 2energy = 0.31701887934635165Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.660092984274442e-05
Newton iteration 3energy = 0.31701887826051933Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.768048486573859e-10
Newton iteration 4energy = 0.3170188782605194Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 3.354660656138913e-16
t = 0.2
Newton iteration 0energy = 0.295494356513846Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.18345925931027088
Newton iteration 1energy = 0.2791640303389012Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.008014651163108922
Newton iteration 2energy = 0.27913186993623224Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.6250158077613474e-05
Newton iteration 3energy = 0.279131869804198Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 6.857117247134904e-11
Newton iteration 4energy = 0.279131869804198Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 3.7828954071306115e-16
t = 0.30000000000000004
Newton iteration 0energy = 0.26692789507920534Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.1334965893056886
Newton iteration 1energy = 0.25820503551107077Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.004188513106581715
Newton iteration 2energy = 0.258196257629024Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.3667568189871645e-06
Newton iteration 3energy = 0.2581962576194896Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.850620349308043e-12
Newton iteration 4energy = 0.25819625761948956Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 3.345785440454107e-16
t = 0.4
Newton iteration 0energy = 0.25204311423689185Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.09236175483868572
Newton iteration 1energy = 0.24783887867786306Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0019843467589302334
Newton iteration 2energy = 0.2478369092241839Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.685217012146076e-07
Newton iteration 3energy = 0.24783690922371493Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 2.34741928405964e-13
Newton iteration 4energy = 0.24783690922371482Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 4 , A u 4 >_{-1}^0.5 = 3.7019931635429417e-16
t = 0.5
Newton iteration 0energy = 0.24499322447001157Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.06164684190451959
Newton iteration 1energy = 0.24311100787346115Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.0008779623506489436
Newton iteration 2energy = 0.24311062240965625Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.8822078759980443e-07
Newton iteration 3energy = 0.24311062240963857Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 8.812738401378623e-15
t = 0.6
Newton iteration 0energy = 0.24187354814794007Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0401663587565137
Newton iteration 1energy = 0.2410718102821203Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.000371536748024391
Newton iteration 2energy = 0.24107174125819195Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.366894131780417e-08
Newton iteration 3energy = 0.24107174125819136Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 5.177430126156407e-16
t = 0.7
Newton iteration 0energy = 0.24055438873764573Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.02578203385924664
Newton iteration 1energy = 0.2402233339155409Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 0.00015331214259509088
Newton iteration 2energy = 0.2402233221629405Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.795968720638723e-09
Newton iteration 3energy = 0.24022332216294048Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.1383921569936535e-16
t = 0.7999999999999999
Newton iteration 0energy = 0.2400118618581594Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.016425332383409572
Newton iteration 1energy = 0.23987730469401894Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 6.283809116183802e-05
Newton iteration 2energy = 0.23987730271968508Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.0154930579076878e-09
Newton iteration 3energy = 0.23987730271968513Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.6394360215478465e-16
t = 0.8999999999999999
Newton iteration 0energy = 0.2397915819363539Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.010464445156609927
Newton iteration 1energy = 0.2397369188905994Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.620620881243329e-05
Newton iteration 2energy = 0.23973691854721507Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.9994608625477876e-10
Newton iteration 3energy = 0.23973691854721502Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.712449072357031e-16
t = 0.9999999999999999
Newton iteration 0energy = 0.23970186094053136Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.006735415826965045
Newton iteration 1energy = 0.2396792028920524Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.1556763084997694e-05
Newton iteration 2energy = 0.23967920282527277Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.099460594642426e-11
Newton iteration 3energy = 0.2396792028252729Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.3603056960577796e-16
t = 1.0999999999999999
Newton iteration 0energy = 0.23966435830447483Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.004452542361403728
Newton iteration 1energy = 0.23965445350959919Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 5.694142646688917e-06
Newton iteration 2energy = 0.23965445349338754Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.772503968582708e-11
Newton iteration 3energy = 0.2396544534933875Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.506057912655227e-16
t = 1.2
Newton iteration 0energy = 0.23964767142737936Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.003097908141817561
Newton iteration 1energy = 0.23964287580252738Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 3.2840536759779764e-06
Newton iteration 2energy = 0.23964287579713486Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 7.767219511989633e-12
Newton iteration 3energy = 0.2396428757971348Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.61045832186731e-16
t = 1.3
Newton iteration 0energy = 0.2396393498836377Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.002326255780763519
Newton iteration 1energy = 0.23963664545354452Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.2181767038734202e-06
Newton iteration 2energy = 0.23963664545108435Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.983596215905177e-12
Newton iteration 3energy = 0.23963664545108426Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.282835092824085e-16
t = 1.4000000000000001
Newton iteration 0energy = 0.2396344792626695Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0019006816819442077
Newton iteration 1energy = 0.2396326736515232Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.687946260850493e-06
Newton iteration 2energy = 0.23963267365009858Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.2982453014628886e-12
Newton iteration 3energy = 0.23963267365009855Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.3568724144101285e-16
t = 1.5000000000000002
Newton iteration 0energy = 0.2396311154900338Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0016629280504841238
Newton iteration 1energy = 0.23962973322336592Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.3855651852122846e-06
Newton iteration 2energy = 0.2396297332224059Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.4586272866410796e-12
Newton iteration 3energy = 0.23962973322240597Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.4042630649306477e-16
t = 1.6000000000000003
Newton iteration 0energy = 0.23962848049094232Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0015187632040499007
Newton iteration 1energy = 0.2396273274180997Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.1887244241389382e-06
Newton iteration 2energy = 0.23962732741739323Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.0010925224110844e-12
Newton iteration 3energy = 0.23962732741739318Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.4047895924105704e-16
t = 1.7000000000000004
Newton iteration 0energy = 0.2396262539867423Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0014193843873494107
Newton iteration 1energy = 0.23962524682330955Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.045632909422474e-06
Newton iteration 2energy = 0.2396252468227629Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 7.300861400170568e-13
Newton iteration 3energy = 0.23962524682276293Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.513659338716933e-16
t = 1.8000000000000005
Newton iteration 0energy = 0.23962429644294805Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0013419591604952733
Newton iteration 1energy = 0.2396233961271811Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 9.33053216472548e-07
Newton iteration 2energy = 0.2396233961267457Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.563795108964307e-13
Newton iteration 3energy = 0.23962339612674574Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.0552927438094465e-16
t = 1.9000000000000006
Newton iteration 0energy = 0.23962254072403263Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0012763226584933272
Newton iteration 1energy = 0.2396217263037931Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 8.399033658894825e-07
Newton iteration 2energy = 0.23962172630344045Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.3696359967523316e-13
Newton iteration 3energy = 0.23962172630344042Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.279057390248914e-16
t = 2.0000000000000004
Newton iteration 0energy = 0.2396209495967649Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0012179715961899576
Newton iteration 1energy = 0.23962020792824323Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 7.604668960240027e-07
Newton iteration 2energy = 0.23962020792795402Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.502994241571864e-13
Newton iteration 3energy = 0.2396202079279541Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.1540976898622703e-16
t = 2.1000000000000005
Newton iteration 0energy = 0.2396194989264497Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0011648358763689138
Newton iteration 1energy = 0.23961882055016812Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 6.914942047541799e-07
Newton iteration 2energy = 0.23961882054992892Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.849344039398136e-13
Newton iteration 3energy = 0.23961882054992895Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.022060791912109e-16
t = 2.2000000000000006
Newton iteration 0energy = 0.2396181709174179Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.001115880689979333
Newton iteration 1energy = 0.23961754835794208Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 6.309368258332487e-07
Newton iteration 2energy = 0.23961754835774315Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.342351024452693e-13
Newton iteration 3energy = 0.2396175483577431Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.124775316997188e-16
t = 2.3000000000000007
Newton iteration 0energy = 0.23961695130021216Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0010705168654718392
Newton iteration 1energy = 0.2396163783254567Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 5.773776868941566e-07
Newton iteration 2energy = 0.23961637832529006Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.9415448405273936e-13
Newton iteration 3energy = 0.23961637832529Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.0997221812292823e-16
t = 2.400000000000001
Newton iteration 0energy = 0.23961582808585116Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0010283540627855545
Newton iteration 1energy = 0.23961529935305362Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 5.297627567863049e-07
Newton iteration 2energy = 0.23961529935291334Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.6199821162823026e-13
Newton iteration 3energy = 0.23961529935291334Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.0039485505058225e-16
t = 2.500000000000001
Newton iteration 0energy = 0.2396147909550374Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0009890970729080236
Newton iteration 1energy = 0.2396143018177979Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 4.872666035318775e-07
Newton iteration 2energy = 0.2396143018176792Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.3598853043606795e-13
Newton iteration 3energy = 0.2396143018176792Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.323276194821206e-16
t = 2.600000000000001
Newton iteration 0energy = 0.23961383091454508Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0009525015914657844
Newton iteration 1energy = 0.23961337730105753Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 4.4921981868305184e-07
Newton iteration 2energy = 0.23961337730095666Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.1473288840449472e-13
Newton iteration 3energy = 0.2396133773009567Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 3 , A u 3 >_{-1}^0.5 = 4.1675889114543924e-16
t = 2.700000000000001
Newton iteration 0energy = 0.2396129400737677Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0009183547757921798
Newton iteration 1energy = 0.23961251839967956Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 4.150665657909605e-07
Newton iteration 2energy = 0.2396125183995935Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.724524553808296e-14
t = 2.800000000000001
Newton iteration 0energy = 0.2396121114812304Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0008864662194052857
Newton iteration 1energy = 0.23961171858168248Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 3.843377621229476e-07
Newton iteration 2energy = 0.23961171858160873Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.279056742710236e-14
t = 2.9000000000000012
Newton iteration 0energy = 0.2396113389953494Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0008566633749443014
Newton iteration 1energy = 0.2396109720692427Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 3.5663295477607706e-07
Newton iteration 2energy = 0.23961097206917908Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 7.079233340119339e-14
t = 3.0000000000000013
Newton iteration 0energy = 0.2396106171776509Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0008287889263646606
Newton iteration 1energy = 0.23961027374068344Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 3.316073784868388e-07
Newton iteration 2energy = 0.2396102737406285Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 6.084230609935835e-14
t = 3.1000000000000014
Newton iteration 0energy = 0.23960994120243975Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.000802699050678995
Newton iteration 1energy = 0.23960961904696895Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 3.0896230208916104e-07
Newton iteration 2energy = 0.2396096190469213Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 5.2419355816600374e-14
t = 3.2000000000000015
Newton iteration 0energy = 0.23960930677943418Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0007782621114844118
Newton iteration 1energy = 0.23960900393990864Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.8843758042950464e-07
Newton iteration 2energy = 0.239609003939867Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 4.53837460691867e-14
t = 3.3000000000000016
Newton iteration 0energy = 0.23960871008705956Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0007553575814402505
Newton iteration 1energy = 0.2396084248101213Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.6980575459079707e-07
Newton iteration 2energy = 0.23960842481008487Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.943490918739899e-14
t = 3.4000000000000017
Newton iteration 0energy = 0.239608147714714Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0007338751008235469
Newton iteration 1energy = 0.23960787843327308Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.528672818809477e-07
Newton iteration 2energy = 0.23960787843324113Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.4395575841844475e-14
t = 3.5000000000000018
Newton iteration 0energy = 0.2396076166126835Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0007137136279143667
Newton iteration 1energy = 0.23960736192339643Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.3744661572988932e-07
Newton iteration 2energy = 0.23960736192336823Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 3.010172589440937e-14
t = 3.600000000000002
Newton iteration 0energy = 0.23960711404861693Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0006947806591329083
Newton iteration 1energy = 0.23960687269229805Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.2338894231226247e-07
Newton iteration 2energy = 0.23960687269227302Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.6431184772178026e-14
t = 3.700000000000002
Newton iteration 0energy = 0.23960663756966022Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0006769915073033805
Newton iteration 1energy = 0.2396064084142239Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 2.1055743624060828e-07
Newton iteration 2energy = 0.2396064084142018Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.3344047294241317e-14
t = 3.800000000000002
Newton iteration 0energy = 0.23960618496947494Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0006602686315022897
Newton iteration 1energy = 0.23960596699506795Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.988309330958618e-07
Newton iteration 2energy = 0.2396059669950483Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 2.0667323036294112e-14
t = 3.900000000000002
Newton iteration 0energy = 0.23960575425948313Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0006445410144724198
Newton iteration 1energy = 0.23960554654551278Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.881019423294816e-07
Newton iteration 2energy = 0.2396055465454951Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.8392439447062177e-14
t = 4.000000000000002
Newton iteration 0energy = 0.239605343643775Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0006297435848349288
Newton iteration 1energy = 0.2396051453575803Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.7827494164571733e-07
Newton iteration 2energy = 0.23960514535756447Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.6359534775821174e-14
t = 4.100000000000001
Newton iteration 0energy = 0.23960495149718894Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.000615816681932275
Newton iteration 1energy = 0.2396047618841374Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.6926490629353927e-07
Newton iteration 2energy = 0.23960476188412316Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.4650955748249466e-14
t = 4.200000000000001
Newton iteration 0energy = 0.23960457634614438Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.000602705561393757
Newton iteration 1energy = 0.23960439472096595Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.6099603651134076e-07
Newton iteration 2energy = 0.23960439472095302Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.316632014806186e-14
t = 4.300000000000001
Newton iteration 0energy = 0.23960421685186575Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005903599395996324
Newton iteration 1energy = 0.23960404259105958Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.534006526259395e-07
Newton iteration 2energy = 0.2396040425910478Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.1891481034967456e-14
t = 4.4
Newton iteration 0energy = 0.2396038717956791Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005787335752184119
Newton iteration 1energy = 0.23960370433085504Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.4641823370105469e-07
Newton iteration 2energy = 0.23960370433084427Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 1.0716509618600795e-14
t = 4.5
Newton iteration 0energy = 0.23960354006611384Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005677838859552442
Newton iteration 1energy = 0.23960337887814379Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.3999457888852992e-07
Newton iteration 2energy = 0.23960337887813393Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 9.745816221516653e-15
t = 4.6
Newton iteration 0energy = 0.2396032206475711Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005574715986072078
Newton iteration 1energy = 0.23960306526144695Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.3408107469814374e-07
Newton iteration 2energy = 0.2396030652614379Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.896636967396773e-15
t = 4.699999999999999
Newton iteration 0energy = 0.23960291261035707Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005477604304859979
Newton iteration 1energy = 0.23960276259066243Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.2863405339167213e-07
Newton iteration 2energy = 0.23960276259065416Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 8.158359737324676e-15
t = 4.799999999999999
Newton iteration 0energy = 0.23960261510190303Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005386168002506355
Newton iteration 1energy = 0.2396024700488202Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.2361423087280676e-07
Newton iteration 2energy = 0.2396024700488126Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 7.475013030957406e-15
t = 4.899999999999999
Newton iteration 0energy = 0.23960232733902043Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 0 , A u 0 >_{-1}^0.5 = 0.0005300095661919717
Newton iteration 1energy = 0.2396021868848042Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 1 , A u 1 >_{-1}^0.5 = 1.1898621270787481e-07
Newton iteration 2energy = 0.2396021868847971Assemble linearization
assemble VOL element 58/58
call pardiso ... done
<A u 2 , A u 2 >_{-1}^0.5 = 6.9022263020530444e-15
t = 4.999999999999998
Minimal energy extension (postscript in unit-2.1.3 )¶
[14]:
from netgen.geom2d import unit_square
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
fes_ho = Discontinuous(H1(mesh, order=10))
fes_lo = H1(mesh, order=1, dirichlet=".*")
fes_lam = Discontinuous(H1(mesh, order=1))
fes = FESpace([fes_ho, fes_lo, fes_lam])
uho, ulo, lam = fes.TrialFunction()
Generate Mesh from spline geometry
Boundary mesh done, np = 40
CalcLocalH: 40 Points 0 Elements 0 Surface Elements
Meshing domain 1 / 1
Surface meshing done
Edgeswapping, topological
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Edgeswapping, metric
Smoothing
Split improve
Combine improve
Smoothing
Update mesh topology
Update clusters
[15]:
a = BilinearForm(fes)
a += Variation(0.5 * grad(uho)*grad(uho)*dx
- 1*uho*dx
+ (uho-ulo)*lam*dx(element_vb=BBND))
gfu = GridFunction(fes)
solvers.Newton(a=a, u=gfu)
Draw(gfu.components[0])
Newton iteration 0
Assemble linearization
assemble VOL element 232/232
call umfpack ... done
err = 0.39392955469293095
Newton iteration 1
Assemble linearization
assemble VOL element 232/232
call umfpack ... done
err = 1.7898870223932246e-15
The minimization problem is solved by the solution of the PDE:
under the constraint
[ ]: