예제 #1
0
 def setUp(self):
     super().setUp()
     domain, geom = mesh.rectilinear([numpy.linspace(0, 1, 9)] * 2)
     ubasis = domain.basis('std', degree=2)
     if self.vector:
         u = ubasis.vector(2).dot(
             function.Argument('dofs', [len(ubasis) * 2]))
     else:
         u = (ubasis[:, numpy.newaxis] *
              function.Argument('dofs', [len(ubasis), 2])).sum(0)
     Geom = geom * [1.1, 1] + u
     self.cons = solver.optimize('dofs',
                                 domain.boundary['left,right'].integral(
                                     (u**2).sum(0), degree=4),
                                 droptol=1e-15)
     self.boolcons = ~numpy.isnan(self.cons)
     strain = .5 * (function.outer(Geom.grad(geom), axis=1).sum(0) -
                    function.eye(2))
     self.energy = domain.integral(
         ((strain**2).sum([0, 1]) + 20 *
          (function.determinant(Geom.grad(geom)) - 1)**2) *
         function.J(geom),
         degree=6)
     self.residual = self.energy.derivative('dofs')
     self.tol = 1e-10
예제 #2
0
 def test_nonlinear_diagonalshift(self):
     nelems = 10
     domain, geom = mesh.rectilinear([nelems, 1])
     geom *= [2 * numpy.pi / nelems, 1]
     ubasis = domain.basis('spline', degree=2).vector(2)
     u = ubasis.dot(function.Argument('dofs', [len(ubasis)]))
     Geom = [.5 * geom[0], geom[1] + function.cos(geom[0])
             ] + u  # compress by 50% and buckle
     cons = solver.minimize('dofs',
                            domain.boundary['left,right'].integral(
                                (u**2).sum(0), degree=4),
                            droptol=1e-15).solve()
     strain = .5 * (function.outer(Geom.grad(geom), axis=1).sum(0) -
                    function.eye(2))
     energy = domain.integral(
         ((strain**2).sum([0, 1]) + 150 *
          (function.determinant(Geom.grad(geom)) - 1)**2) *
         function.J(geom),
         degree=6)
     nshift = 0
     for iiter, (lhs, info) in enumerate(
             solver.minimize('dofs', energy, constrain=cons)):
         self.assertLess(iiter, 38)
         if info.shift:
             nshift += 1
         if info.resnorm < self.tol:
             break
     self.assertEqual(nshift, 9)
예제 #3
0
 def setUp(self):
   super().setUp()
   domain, geom = mesh.rectilinear([numpy.linspace(0,1,9)] * 2)
   ubasis = domain.basis('std', degree=2).vector(2)
   u = ubasis.dot(function.Argument('dofs', [len(ubasis)]))
   Geom = geom * [1.1, 1] + u
   self.cons = solver.minimize('dofs', domain.boundary['left,right'].integral((u**2).sum(0), degree=4), droptol=1e-15).solve()
   self.boolcons = ~numpy.isnan(self.cons)
   strain = .5 * (function.outer(Geom.grad(geom), axis=1).sum(0) - function.eye(2))
   self.energy = domain.integral(((strain**2).sum([0,1]) + 20*(function.determinant(Geom.grad(geom))-1)**2)*function.J(geom), degree=6)
   self.residual = self.energy.derivative('dofs')
   self.tol = 1e-10
예제 #4
0
 def test_nonlinear_diagonalshift(self):
   nelems = 10
   domain, geom = mesh.rectilinear([nelems,1])
   geom *= [2*numpy.pi/nelems, 1]
   ubasis = domain.basis('spline', degree=2).vector(2)
   u = ubasis.dot(function.Argument('dofs', [len(ubasis)]))
   Geom = [.5 * geom[0], geom[1] + function.cos(geom[0])] + u # compress by 50% and buckle
   cons = solver.minimize('dofs', domain.boundary['left,right'].integral((u**2).sum(0), degree=4), droptol=1e-15).solve()
   strain = .5 * (function.outer(Geom.grad(geom), axis=1).sum(0) - function.eye(2))
   energy = domain.integral(((strain**2).sum([0,1]) + 150*(function.determinant(Geom.grad(geom))-1)**2)*function.J(geom), degree=6)
   nshift = 0
   for iiter, (lhs, info) in enumerate(solver.minimize('dofs', energy, constrain=cons)):
     self.assertLess(iiter, 90)
     if info.shift:
       nshift += 1
     if info.resnorm < self.tol:
       break
   self.assertEqual(nshift, 60)
예제 #5
0
    with _plot('p', mesh=mesh, **kwargs) as (fig, ax):
        im = ax.tripcolor(tri, pvals, shading='gouraud')
        _colorbar(fig, im, **kwargs)


def deformation(case, mu, lhs, stress='xy', name='solution', **kwargs):
    disp = case.bases['u'].obj.dot(lhs)
    refgeom = case.geometry(mu)
    geom = refgeom + disp

    E = mu['ymod1']
    NU = mu['prat']
    MU = E / (1 + NU)
    LAMBDA = E * NU / (1 + NU) / (1 - 2 * NU)
    stressfunc = -MU * disp.symgrad(refgeom) + LAMBDA * disp.div(
        refgeom) * fn.eye(disp.shape[0])

    if geom.shape == (2, ):
        stressfunc = stressfunc[tuple('xyz'.index(c) for c in stress)]
        mesh, stressdata = case.domain.elem_eval([geom, 1],
                                                 separate=True,
                                                 ischeme='bezier3')
        with _plot('u', name=name, **kwargs) as plt:
            plt.mesh(mesh, stressdata)
            _colorbar(plt, **kwargs)

    elif geom.shape == (3, ):
        nutils.plot.writevtu(name,
                             case.domain,
                             geom,
                             pointdata={
예제 #6
0
def main(fname: str, degree: int, Δuload: unit['mm'], nsteps: int, E: unit['GPa'], nu: float, σyield: unit['MPa'], hardening: unit['GPa'], referencedata: typing.Optional[str], testing: bool):

  '''
  Plastic deformation of a perforated strip.

  .. arguments::

     fname [strip.msh]
       Mesh file with units in [mm]

     degree [1]
       Finite element interpolation order

     Δuload [5μm]
       Load boundary displacement steps

     nsteps [11]
       Number of load steps

     E [70GPa]
       Young's modulus

     nu [0.2]
       Poisson ratio

     σyield [243MPa]
       Yield strength

     hardening [2.25GPa]
       Hardening parameter

     referencedata [zienkiewicz.csv]
       Reference data file name

     testing [False]
       Use a 1 element mesh for testing

  .. presets::

     testing
       nsteps=30
       referencedata=
       testing=True
  '''

  # We commence with reading the mesh from the specified GMSH file, or, alternatively,
  # we continue with a single-element mesh for testing purposes.
  domain, geom = mesh.gmsh(pathlib.Path(__file__).parent/fname)
  Wdomain, Hdomain = domain.boundary['load'].integrate([function.J(geom),geom[1]*function.J(geom)], degree=1)
  Hdomain /= Wdomain

  if testing:
    domain, geom = mesh.rectilinear([numpy.linspace(0,Wdomain/2,2),numpy.linspace(0,Hdomain,2)])
    domain = domain.withboundary(hsymmetry='left', vsymmetry='bottom', load='top')

  # We next initiate the point set in which the constitutive bahvior will be evaluated.
  # Note that this is also the point set in which the history variable will be stored.
  gauss  = domain.sample('gauss', 2)

  # Elasto-plastic formulation
  # ==========================
  # The weak formulation is constructed using a `Namespace`, which is initialized and
  # populated with the necessary problem parameters and coordinate system. Note that the
  # coefficients `mu` and `λ` have been defined such that 'C_{ijkl}' is the fourth-order
  # **plane stress** elasticity tensor.

  ns = function.Namespace(fallback_length=domain.ndims)

  ns.x      = geom
  ns.Δuload = Δuload
  ns.mu     = E/(1-nu**2)*((1-nu)/2)
  ns.λ      = E/(1-nu**2)*nu
  ns.delta  = function.eye(domain.ndims)
  ns.C_ijkl = 'mu ( delta_ik delta_jl + delta_il delta_jk ) + λ delta_ij delta_kl'

  # We make use of a Lagrange finite element basis of arbitrary `degree`. Since we
  # approximate the displacement field, the basis functions are vector-valued. Both
  # the displacement field at the current load step `u` and displacement field of the
  # previous load step `u0` are approximated using this basis:

  ns.basis = domain.basis('std',degree=degree).vector(domain.ndims)

  ns.u0_i = 'basis_ni ?lhs0_n'
  ns.u_i  = 'basis_ni ?lhs_n'

  # This simulation is based on a standard elasto-plasticity model. In this
  # model the *total strain*, ε_kl = ½ (∂u_k/∂x_l + ∂u_l/∂x_k)', is comprised of an
  # elastic and a plastic part:
  #
  #   ε_kl = εe_kl + εp_kl
  #
  # The stress is related to the *elastic strain*, εe_kl, through Hooke's law:
  #
  #   σ_ij = C_ijkl εe_kl = C_ijkl (ε_kl - εp_kl)

  ns.ε_kl   = '(u_k,l + u_l,k) / 2'
  ns.ε0_kl  = '(u0_k,l + u0_l,k) / 2'
  ns.gbasis = gauss.basis()
  ns.εp0_ij = 'gbasis_n ?εp0_nij'
  ns.κ0     = 'gbasis_n ?κ0_n'

  ns.εp    = PlasticStrain(ns.ε, ns.ε0, ns.εp0, ns.κ0, E, nu, σyield, hardening)
  ns.εe_ij = 'ε_ij - εp_ij'
  ns.σ_ij  = 'C_ijkl εe_kl'

  # Note that the plasticity model is implemented through the user-defined function `PlasticStrain`,
  # which implements the actual yielding model including a standard return mapping algorithm. This
  # function is discussed in detail below.
  #
  # The components of the residual vector are then defined as:
  #
  #   r_n = ∫_Ω (∂N_ni/∂x_j) σ_ij dΩ

  res = domain.integral('basis_ni,j σ_ij d:x' @ ns, degree=2)

  # The problem formulation is completed by supplementing prescribed displacement boundary conditions
  # (for a load step), which are computed in the standard manner:

  sqr  = domain.boundary['hsymmetry,vsymmetry'].integral('(u_k n_k)^2 d:x' @ ns, degree=2)
  sqr += domain.boundary['load'].integral('(u_k n_k - Δuload)^2 d:x' @ ns, degree=2)
  cons = solver.optimize('lhs', sqr, droptol=1e-15)

  # Incremental-iterative solution procedure
  # ========================================
  # We initialize the solution vector for the first load step `lhs` and solution vector
  # of the previous load step `lhs0`. Note that, in order to construct a predictor step
  # for the first load step, we define the previous load step state as the solution
  # vector corresponding to a negative elastic loading step.

  lhs0 = -solver.solve_linear('lhs', domain.integral('basis_ni,j C_ijkl ε_kl d:x' @ ns, degree=2), constrain=cons)
  lhs  = numpy.zeros_like(lhs0)
  εp0  = numpy.zeros((gauss.npoints,)+ns.εp0.shape)
  κ0   = numpy.zeros((gauss.npoints,)+ns.κ0.shape)

  # To store the force-dispalcement data we initialize an empty data array with the
  # inital state solution substituted in the first row.
  fddata      = numpy.empty(shape=(nsteps+1,2))
  fddata[:]   = numpy.nan
  fddata[0,:] = 0

  # Load step incrementation
  # ------------------------
  with treelog.iter.fraction('step', range(nsteps)) as counter:
    for step in counter:

      # The solution of the previous load step is set to `lhs0`, and the Newton solution
      # procedure is initialized by extrapolation of the state vector:
      lhs_init = lhs + (lhs-lhs0)
      lhs0     = lhs

      # The non-linear system of equations is solved for `lhs` using Newton iterations,
      # where the `step` variable is used to scale the incremental constraints.
      lhs = solver.newton(target='lhs', residual=res, constrain=cons*step, lhs0=lhs_init, arguments={'lhs0':lhs0,'εp0':εp0,'κ0':κ0}).solve(tol=1e-6)

      # The computed solution is post-processed in the form of a loading curve - which
      # plots the normalized mean stress versus the maximum 'ε_{yy}' strain
      # component - and a contour plot showing the 'σ_{yy}' stress component on a
      # deformed mesh. Note that since the stresses are defined in the integration points
      # only, a post-processing step is involved that transfers the stress information to
      # the nodal points.
      εyymax = gauss.eval(ns.ε[1,1], arguments=dict(lhs=lhs)).max()

      basis = domain.basis('std', degree=1)
      bw, b = domain.integrate([basis * ns.σ[1,1] * function.J(geom), basis * function.J(geom)], degree=2, arguments=dict(lhs=lhs,lhs0=lhs0,εp0=εp0,κ0=κ0))
      σyy = basis.dot(bw / b)

      uyload, σyyload = domain.boundary['load'].integrate(['u_1 d:x'@ns,σyy * function.J(geom)], degree=2, arguments=dict(lhs=lhs,εp0=εp0,κ0=κ0))
      uyload /= Wdomain
      σyyload /= Wdomain
      fddata[step,0] = (E*εyymax)/σyield
      fddata[step,1] = (σyyload*2)/σyield

      with export.mplfigure('forcedisp.png') as fig:
        ax = fig.add_subplot(111, xlabel=r'${E \cdot {\rm max}(\varepsilon_{yy})}/{\sigma_{\rm yield}}$', ylabel=r'${\sigma_{\rm mean}}/{\sigma_{\rm yield}}$')
        if referencedata:
          data = numpy.genfromtxt(pathlib.Path(__file__).parent/referencedata, delimiter=',', skip_header=1)
          ax.plot(data[:,0], data[:,1], 'r:', label='Reference')
          ax.legend()
        ax.plot(fddata[:,0], fddata[:,1], 'o-', label='Nutils')
        ax.grid()

      bezier = domain.sample('bezier', 3)
      points, uvals, σyyvals = bezier.eval(['(x_i + 25 u_i)' @ ns, ns.u, σyy], arguments=dict(lhs=lhs))
      with export.mplfigure('stress.png') as fig:
        ax = fig.add_subplot(111, aspect='equal', xlabel=r'$x$ [mm]', ylabel=r'$y$ [mm]')
        im = ax.tripcolor(points[:,0]/unit('mm'), points[:,1]/unit('mm'), bezier.tri, σyyvals/unit('MPa'), shading='gouraud', cmap='jet')
        ax.add_collection(collections.LineCollection(points.take(bezier.hull, axis=0), colors='k', linewidths=.1))
        ax.autoscale(enable=True, axis='both', tight=True)
        cb = fig.colorbar(im)
        im.set_clim(0, 1.2*σyield)
        cb.set_label(r'$σ_{yy}$ [MPa]')

      # Load step convergence
      # ---------------------
      # At the end of the loading step, the plastic strain state and history parameter are updated,
      # where use if made of the strain hardening relation for the history variable:
      #
      #   Δκ = √(Δεp_ij Δεp_ij)

      Δεp = ns.εp-ns.εp0
      Δκ  = function.sqrt((Δεp*Δεp).sum((0,1)))
      κ0  = gauss.eval(ns.κ0+Δκ, arguments=dict(lhs0=lhs0,lhs=lhs,εp0=εp0,κ0=κ0))
      εp0 = gauss.eval(ns.εp, arguments=dict(lhs0=lhs0,lhs=lhs,εp0=εp0,κ0=κ0))
예제 #7
0
def main(
    nelems: 'number of elements' = 12,
    viscosity: 'fluid viscosity' = 1e-2,
    density: 'fluid density' = 1,
    tol: 'solver tolerance' = 1e-12,
    rotation: 'cylinder rotation speed' = 0,
    timestep: 'time step' = 1/24,
    maxradius: 'approximate domain size' = 25,
    tmax: 'end time' = numpy.inf,
    withplots: 'create plots' = True,
  ):

  uinf = numpy.array([ 1, 0 ])
  log.user( 'reynolds number:', density/viscosity )

  # construct mesh
  rscale = numpy.pi / nelems
  melems = numpy.ceil( numpy.log(2*maxradius) / rscale ).astype( int )
  log.info( 'creating {}x{} mesh, outer radius {:.2f}'.format( melems, 2*nelems, .5*numpy.exp(rscale*melems) ) )
  domain, gridgeom = mesh.rectilinear( [range(melems+1),numpy.linspace(0,2*numpy.pi,2*nelems+1)], periodic=(1,) )
  rho, phi = gridgeom
  phi += 1e-3 # tiny nudge (0.057 deg) to break element symmetry
  cylvelo = -.5 * rotation * function.trignormal(phi)
  radius = .5 * function.exp( rscale * rho )
  geom = radius * function.trigtangent(phi)

  # prepare bases
  J = geom.grad( gridgeom )
  detJ = function.determinant( J )
  vnbasis, vtbasis, pbasis = function.chain([ # compatible spaces using piola transformation
    domain.basis( 'spline', degree=(3,2), removedofs=((0,),None) )[:,_] * J[:,0] / detJ,
    domain.basis( 'spline', degree=(2,3) )[:,_] * J[:,1] / detJ,
    domain.basis( 'spline', degree=2 ) / detJ,
  ])
  vbasis = vnbasis + vtbasis
  stressbasis = (2*viscosity) * vbasis.symgrad(geom) - pbasis[:,_,_] * function.eye( domain.ndims )

  # prepare matrices
  A = function.outer( vbasis.grad(geom), stressbasis ).sum([2,3]) + function.outer( pbasis, vbasis.div(geom) )
  Ao = function.outer( vbasis, density * ( vbasis.grad(geom) * uinf ).sum(-1) ).sum(-1)
  At = (density/timestep) * function.outer( vbasis ).sum(-1)
  Ad = function.outer( vbasis.div(geom) )
  stokesmat, uniconvec, inertmat, divmat = domain.integrate( [ A, Ao, At, Ad ], geometry=geom, ischeme='gauss9' )

  # interior boundary condition (weak imposition of shear verlocity component)
  inner = domain.boundary['left']
  h = inner.elem_eval( 1, geometry=geom, ischeme='gauss9', asfunction=True )
  nietzsche = (2*viscosity) * ( (7.5/h) * vbasis - vbasis.symgrad(geom).dotnorm(geom) )
  B = function.outer( nietzsche, vbasis ).sum(-1)
  b = ( nietzsche * cylvelo ).sum(-1)
  bcondmat, rhs = inner.integrate( [ B, b ], geometry=geom, ischeme='gauss9' )
  stokesmat += bcondmat

  # exterior boundary condition (full velocity vector imposed at inflow)
  inflow = domain.boundary['right'].select( -( uinf * geom.normal() ).sum(-1), ischeme='gauss1' )
  cons = inflow.project( uinf, onto=vbasis, geometry=geom, ischeme='gauss9', tol=1e-12 )

  # initial condition (stationary oseen flow)
  lhs = (stokesmat+uniconvec).solve( rhs, constrain=cons, tol=0 )

  # prepare plotting
  if not withplots:
    makeplots = lambda *args: None
  else:
    makeplots = MakePlots( domain, geom, video=withplots=='video', timestep=timestep )
    mesh_ = domain.elem_eval( geom, ischeme='bezier5', separate=True )
    inflow_ = inflow.elem_eval( geom, ischeme='bezier5', separate=True )
    with plot.PyPlot( 'mesh', ndigits=0 ) as plt:
      plt.mesh( mesh_ )
      plt.rectangle( makeplots.bbox[:,0], *( makeplots.bbox[:,1] - makeplots.bbox[:,0] ), ec='green' )
      plt.segments( inflow_, colors='red' )

  # start time stepping
  stokesmat += inertmat
  precon = (stokesmat+uniconvec).getprecon( 'splu', constrain=cons )
  for istep in log.count( 'timestep', start=1 ):
    log.info( 'velocity divergence:', divmat.matvec(lhs).dot(lhs) )
    convection = density * ( vbasis.grad(geom) * vbasis.dot(lhs) ).sum(-1)
    matrix = stokesmat + domain.integrate( function.outer( vbasis, convection ).sum(-1), ischeme='gauss9', geometry=geom )
    lhs = matrix.solve( rhs + inertmat.matvec(lhs), lhs0=lhs, constrain=cons, tol=tol, restart=9999, precon=precon )
    makeplots( vbasis.dot(lhs), pbasis.dot(lhs), istep*timestep*rotation )
    if istep*timestep > tmax:
      break

  return lhs
예제 #8
0
def main(
    nelems: 'number of elements' = 12,
    viscosity: 'fluid viscosity' = 1e-3,
    density: 'fluid density' = 1,
    degree: 'polynomial degree' = 2,
    warp: 'warp domain (downward bend)' = False,
    tol: 'solver tolerance' = 1e-5,
    maxiter: 'maximum number if iterations, 0 for unlimited' = 0,
    withplots: 'create plots' = True,
  ):

  Re = density / viscosity # based on unit length and velocity
  log.user( 'reynolds number: {:.1f}'.format(Re) )

  # construct mesh
  verts = numpy.linspace( 0, 1, nelems+1 )
  domain, geom = mesh.rectilinear( [verts,verts] )

  # construct bases
  vxbasis, vybasis, pbasis, lbasis = function.chain([
    domain.basis( 'spline', degree=(degree+1,degree), removedofs=((0,-1),None) ),
    domain.basis( 'spline', degree=(degree,degree+1), removedofs=(None,(0,-1)) ),
    domain.basis( 'spline', degree=degree ),
    [1], # lagrange multiplier
  ])
  if not warp:
    vbasis = function.stack( [ vxbasis, vybasis ], axis=1 )
  else:
    gridgeom = geom
    xi, eta = gridgeom
    geom = (eta+2) * function.rotmat(xi*.4)[:,1] - (0,2) # slight downward bend
    J = geom.grad( gridgeom )
    detJ = function.determinant( J )
    vbasis = ( vxbasis[:,_] * J[:,0] + vybasis[:,_] * J[:,1] ) / detJ # piola transform
    pbasis /= detJ
  stressbasis = (2*viscosity) * vbasis.symgrad(geom) - (pbasis)[:,_,_] * function.eye( domain.ndims )

  # construct matrices
  A = function.outer( vbasis.grad(geom), stressbasis ).sum([2,3]) \
    + function.outer( pbasis, vbasis.div(geom)+lbasis ) \
    + function.outer( lbasis, pbasis )
  Ad = function.outer( vbasis.div(geom) )
  stokesmat, divmat = domain.integrate( [ A, Ad ], geometry=geom, ischeme='gauss9' )

  # define boundary conditions
  normal = geom.normal()
  utop = function.asarray([ normal[1], -normal[0] ])
  h = domain.boundary.elem_eval( 1, geometry=geom, ischeme='gauss9', asfunction=True )
  nietzsche = (2*viscosity) * ( ((degree+1)*2.5/h) * vbasis - vbasis.symgrad(geom).dotnorm(geom) )
  stokesmat += domain.boundary.integrate( function.outer( nietzsche, vbasis ).sum(-1), geometry=geom, ischeme='gauss9' )
  rhs = domain.boundary['top'].integrate( ( nietzsche * utop ).sum(-1), geometry=geom, ischeme='gauss9' )

  # prepare plotting
  makeplots = MakePlots( domain, geom ) if withplots else lambda *args: None

  # start picard iterations
  lhs = stokesmat.solve( rhs, tol=tol, solver='cg', precon='spilu' )
  for iiter in log.count( 'picard' ):
    log.info( 'velocity divergence:', divmat.matvec(lhs).dot(lhs) )
    makeplots( vbasis.dot(lhs), pbasis.dot(lhs) )
    ugradu = ( vbasis.grad(geom) * vbasis.dot(lhs) ).sum(-1)
    convection = density * function.outer( vbasis, ugradu ).sum(-1)
    matrix = stokesmat + domain.integrate( convection, ischeme='gauss9', geometry=geom )
    lhs, info = matrix.solve( rhs, lhs0=lhs, tol=tol, info=True, precon='spilu', restart=999 )
    if iiter == maxiter-1 or info.niter == 0:
      break

  return rhs, lhs