def newton(f, x0, maxiter, restol, lintol):
    '''find x such that f(x) = 0, starting at x0'''

    for i in log.range('newton', maxiter + 1):
        if i > 0:
            # solve system linearised around `x`, solution in `direction`
            J = f(x, jacobian=True)
            direction = -J.solve(residual, tol=lintol)
            xprev = x
            residual_norm_prev = residual_norm
            # line search, stop if there is convergence w.r.t. iteration `i-1`
            for j in range(4):
                x = xprev + (0.5)**j * direction
                residual = f(x)
                residual_norm = numpy.linalg.norm(residual, 2)
                if residual_norm < residual_norm_prev:
                    break
                # else: no convergence, repeat and reduce step size by one half
            else:
                log.warning('divergence')
        else:
            # before first iteration, compute initial residual
            x = x0
            residual = f(x)
            residual_norm = numpy.linalg.norm(residual, 2)
        log.info('residual norm: {:.5e}'.format(residual_norm))
        if residual_norm < restol:
            break
    else:
        raise ValueError(
            'newton did not converge in {} iterations'.format(maxiter))
    return x
Ejemplo n.º 2
0
    def __init__(self, g, order, f=None, **kwargs):
        if order < 6:
            log.warning('''
                    Warning, for this method a Gauss-scheme of at least
                    order 6 is recommended.
                ''')
        super().__init__(g, order, **kwargs)

        if isinstance(f, go.TensorGridObject):

            if len(f) != 2:
                raise NotImplementedError

            if not (f.knotvector <= g.knotvector).all():
                log.warning('''
                        Warning, using a control mapping whose knotvector is
                        not a subset of the target GridObject knotvector
                        or that does not have the same periodicity properties
                        is not recommended, since Gaussian quadrature is
                        ill-defined in this case.
                    ''')

            self._f = f.toscipy()

        else:
            '''
                So far, we only allow for instantiations of the class
                ``go.TensorGridObject`` as control mapping.
                XXX: allow for other means of instantiation, for instance
                via a list / dictionary implementing functions for the
                control mapping and all its relevant derivatives.
            '''
            raise NotImplementedError

        self._tablulate_control_mapping()
Ejemplo n.º 3
0
    def from_parent(cls, parent, *key):
        assert isinstance(parent, cls)
        if len(key) == 1:
            side, = key
            d = parent.__dict__.copy()
            d['domain'] = d['domain'].boundary[side]

            # flatten appropriate axis
            bnames = parent.domain._bnames
            sidetovalue = dict(zip(bnames, range(len(bnames))))
            deleteindex = sidetovalue[side] // 2
            deleteaxis = d['axesnames'][deleteindex]
            d['knotvector'] = np.delete(d['knotvector'], deleteindex)
            geom = d['geom']
            d[ 'geom' ] = function.stack( [ geom[i] for i in range( len( parent ) ) \
                    if not i == deleteindex ] )
            d['axesnames'] = tuple(
                filter(lambda x: x != deleteaxis, d['axesnames']))
            if len(d['knotvector']) == 0:
                log.warning(
                    'The GridObject has not been implemented with an empty knotvector yet'
                )
                raise NotImplementedError

            ret = cls(**d)
            ret.x = parent[side]
            for name in ret.domain._bnames:
                ret._cons[ret.index.boundary(name).flatten] = ret[name]
            return ret
        ret = parent
        for side in key:
            ret = ret.from_parent(ret, side)
        return ret
Ejemplo n.º 4
0
def smoothen_vector(vec0,
                    dt,
                    method='finite_difference',
                    stop={'T': 0.01},
                    miniter=10):
    """
        Smoothen vec[1: -1] using `method' until stopping criterion has been reached,
        while vec[0] and vec[-1] are held fixed.

        Parameters
        ----------
        vec0 : to-be-smoothened vector
        dt : (initial) timestep, is reduced if the convergence criterion is reached
             before #iterations >= ``miniter``
        method : smoothing method
        stop = {
                'T': (finite difference) if t > T, terminate
                'maxiter': (finite difference) if i > maxiter, terminate
                'vec': if not stop[ 'vec' ]( vec ) terminate
               }
        miniter : minimum number of iterations
    """

    t, i = 0, 0
    vec = vec0.copy()
    d = ChainMap(stop, {'T': np.inf, 'maxiter': 100, 'vec': lambda x: True})
    if method == 'finite_difference':
        N = len(vec)
        dx = 1 / (N - 1)
        fac = dt / (dx**2)
        A = sparse.diags([(-fac * np.ones(N - 2)).tolist() + [0], [1] +
                          ((1 + 2 * fac) * np.ones(N - 2)).tolist() + [1],
                          [0] + (-fac * np.ones(N - 2)).tolist()],
                         [-1, 0, 1]).tocsc()
        A = sparse.linalg.splu(A)
        while True:
            if not all([t < d['T'], i < d['maxiter'], d['vec'](vec)]):
                if i <= miniter:  # timestep too big
                    log.info('Initial timestep too big, reducing to {}'.format(
                        dt / 10))
                    return smoothen_vector(vec0,
                                           dt / 10,
                                           method=method,
                                           stop=stop)
                break
            vec = A.solve(vec)
            t += dt
            i += 1
    else:
        raise "Unknown method '{}'".format(method)

    if d['vec'](vec):
        log.warning('Failed to reach the termination criterion')
    else:
        log.info('Criterion reached at t={} in {} iterations'.format(t, i))

    return vec
Ejemplo n.º 5
0
    def __init__(self, g, order, absc=None):

        # any problems with ``g`` will be caught by the parent __init__
        super().__init__(g, order)

        self._dindices = self._g.dofindices

        # The Jacobian-Determinant is evaluated in
        # the tensor-product points of xi and eta
        # as a constraint to warrant bijectivity of the result
        if isinstance(absc, int):
            # if an integer is passed to absc, we construct
            # abscissae that correspond to and ``absc``-order
            # Gauss quadrature scheme over all elements
            absc = fsol.make_quadrature(g, absc)[0]

        self._absc = absc

        if absc is not None:
            assert all(
                all( [aux.isincreasing(x_), x_[0] >= 0, x_[-1] <= 1] )
                                                    for x_ in absc ), \
                'Invalid constraint abscissae received.'

            _bsplev = lambda i, **kwargs: sparse.csr_matrix(
                self._splev(i, self._absc, **kwargs).ravel()[:, None])

            structure = sparse.hstack([_bsplev(i) for i in range(self._N)])
            self._w_xi = sparse.hstack(
                [_bsplev(i, dx=1) for i in range(self._N)])
            self._w_eta = sparse.hstack(
                [_bsplev(i, dy=1) for i in range(self._N)])

            # sparsity structure of the constraint jacobian
            self.jacobianstructure = \
                sparse.hstack( [ structure, structure ] ).\
                tolil()[:, self._dindices ].nonzero()

        log.info('Constraint jacobian sparsity pattern stored.')

        if not (self.constraint(self._g[self._dindices]) >= 0).all():
            log.warning('''
                    Warning, the initial guess corresponding to the passed
                    GridObject is not feasible. Unless another initial guess
                    is specified, the optimizer will most likely fail.
                ''')
Ejemplo n.º 6
0
def main(
    nelems: 'number of elements' = 20,
    epsilon: 'epsilon, 0 for automatic (based on nelems)' = 0,
    timestep: 'time step' = .01,
    maxtime: 'end time' = 1.,
    theta: 'contact angle (degrees)' = 90,
    init: 'initial condition (random/bubbles)' = 'random',
    figures: 'create figures' = True,
):

    mineps = 1. / nelems
    if not epsilon:
        log.info('setting epsilon={}'.format(mineps))
        epsilon = mineps
    elif epsilon < mineps:
        log.warning('epsilon under crititical threshold: {} < {}'.format(
            epsilon, mineps))

    # create namespace
    ns = function.Namespace()
    ns.epsilon = epsilon
    ns.ewall = .5 * numpy.cos(theta * numpy.pi / 180)

    # construct mesh
    xnodes = ynodes = numpy.linspace(0, 1, nelems + 1)
    domain, ns.x = mesh.rectilinear([xnodes, ynodes])

    # prepare bases
    ns.cbasis, ns.mubasis = function.chain(
        [domain.basis('spline', degree=2),
         domain.basis('spline', degree=2)])

    # polulate namespace
    ns.c = 'cbasis_n ?lhs_n'
    ns.c0 = 'cbasis_n ?lhs0_n'
    ns.mu = 'mubasis_n ?lhs_n'
    ns.f = '(6 c0 - 2 c0^3 - 4 c) / epsilon^2'

    # construct initial condition
    if init == 'random':
        numpy.random.seed(0)
        lhs0 = numpy.random.normal(0, .5, ns.cbasis.shape)
    elif init == 'bubbles':
        R1 = .25
        R2 = numpy.sqrt(.5) * R1  # area2 = .5 * area1
        ns.cbubble1 = function.tanh(
            (R1 - function.norm2(ns.x -
                                 (.5 + R2 / numpy.sqrt(2) + .8 * ns.epsilon)))
            / ns.epsilon)
        ns.cbubble2 = function.tanh(
            (R2 - function.norm2(ns.x -
                                 (.5 - R1 / numpy.sqrt(2) - .8 * ns.epsilon)))
            / ns.epsilon)
        sqr = domain.integral('(c - cbubble1 - cbubble2 - 1)^2 + mu^2' @ ns,
                              geometry=ns.x,
                              degree=4)
        lhs0 = solver.optimize('lhs', sqr)
    else:
        raise Exception('unknown init %r' % init)

    # construct residual
    res = domain.integral(
        'epsilon^2 cbasis_n,k mu_,k + mubasis_n (mu + f) - mubasis_n,k c_,k'
        @ ns,
        geometry=ns.x,
        degree=4)
    res -= domain.boundary.integral('mubasis_n ewall' @ ns,
                                    geometry=ns.x,
                                    degree=4)
    inertia = domain.integral('cbasis_n c' @ ns, geometry=ns.x, degree=4)

    # solve time dependent problem
    nsteps = numeric.round(maxtime / timestep)
    makeplots = MakePlots(domain, nsteps,
                          ns) if figures else lambda *args: None
    for istep, lhs in log.enumerate(
            'timestep',
            solver.impliciteuler('lhs',
                                 target0='lhs0',
                                 residual=res,
                                 inertia=inertia,
                                 timestep=timestep,
                                 lhs0=lhs0)):
        makeplots(lhs)
        if istep == nsteps:
            break

    return lhs0, lhs
Ejemplo n.º 7
0
    def _tablulate_control_mapping(self):
        '''
            Tabulate the control mapping and all its required derivatives
            in the quadrature points.
            XXX: see if we can make this look prettier
                (possibly using some symbolic approach).
        '''

        # by the time this is called, _f and _quad have been set
        # in the __init__ method
        _f = self._f
        _q = self._quad

        # all required first and second derivatives
        X_xi = _f(*_q, dx=1)
        X_xi_eta = _f(*_q, dx=1, dy=1)
        X_xi_xi = _f(*_q, dx=2)
        X_eta = _f(*_q, dy=1)
        X_eta_eta = _f(*_q, dy=2)

        _0 = (Ellipsis, 0)
        _1 = (Ellipsis, 1)

        # Jacobian determinant and its derivatives
        jacdet = X_xi[_0] * X_eta[_1] - X_eta[_0] * X_xi[_1]
        jacsqrt = jacdet**2
        jacdet_xi = X_xi_xi[_0] * X_eta[_1] + X_xi[_0] * X_xi_eta[_1] \
            - X_xi_eta[_0] * X_xi[_1] - X_eta[_0] * X_xi_xi[_1]
        jacdet_eta = X_xi_eta[_0] * X_eta[_1] + X_xi[_0] * X_eta_eta[_1] \
            - X_eta_eta[_0] * X_xi[_1] - X_eta[_0] * X_xi_eta[_1]

        if not (jacdet > 0).all():
            log.warning('''
                    Warning, the control mapping is defective on at least one
                    quadrature point. Proceeding with this control mapping may
                    lead to failure of convergence and / or a defective mapping.
                ''')

        # metric tensor of the control mapping divided by the Jacobian determinant
        self._G11 = (X_xi**2).sum(-1) / jacdet
        self._G12 = (X_xi * X_eta).sum(-1) / jacdet
        self._G22 = (X_eta**2).sum(-1) / jacdet

        # derivatives of _G11
        self._G11_xi = \
            2 * ( X_xi * X_xi_xi ).sum(-1) / jacdet \
            - jacdet_xi * (X_xi ** 2).sum(-1) / jacsqrt
        self._G11_eta = \
            2 * ( X_xi * X_xi_eta ).sum(-1) / jacdet \
            - jacdet_eta * (X_xi ** 2).sum(-1) / jacsqrt

        # derivatives of _G12
        self._G12_xi = \
            (X_xi_xi * X_eta + X_xi * X_xi_eta).sum(-1) / jacdet  \
            - jacdet_xi * (X_xi * X_eta).sum(-1) / jacsqrt
        self._G12_eta = \
            (X_xi_eta * X_eta + X_xi * X_eta_eta).sum(-1) / jacdet \
            - jacdet_eta * (X_xi * X_eta).sum(-1) / jacsqrt

        # derivatives of _G22
        self._G22_xi = \
            2 * ( X_eta * X_xi_eta ).sum(-1) / jacdet \
            - jacdet_xi * (X_eta ** 2).sum(-1) / jacsqrt
        self._G22_eta = \
            2 * ( X_eta * X_eta_eta ).sum(-1) / jacdet \
            - jacdet_eta * (X_eta ** 2).sum(-1) / jacsqrt

        log.info('The control mapping has been tabulated.')
Ejemplo n.º 8
0
 def inverse(self):
     ''' Return exact inverse of self '''
     log.warning(
         'Computing the exact inverse is slow, use approximate_inverse instead'
     )
     return self.__class__(pynverse.inversefunc(self._f, domain=[0, 1]))
Ejemplo n.º 9
0
def InterpolatedUnivariateSpline(pc, k=3, c0_thresh=None, **scipyargs):
    if isinstance(pc, np.ndarray):
        pc = PointCloud(pc)
    assert isinstance(pc, PointCloud)

    if k > len(pc) + 1:
        log.warning(
            'Warning, the number of points is too small for a {}-th order interpolations, \
                it will be clipped.'.format(k))
        k = min(k, len(pc) + 1)

    if c0_thresh is None:
        periodic = pc.periodic

        def f(x):
            verts, points = pc.verts, pc.points
            if periodic:  # repeat points to acquire (nearly) continuous derivative at x = 0
                verts = np.concatenate([verts - 1, verts, verts + 1])
                points = np.vstack([points] * 3)
            splines = (interpolate.InterpolatedUnivariateSpline(verts,
                                                                i,
                                                                k=k,
                                                                **scipyargs)
                       for i in points.T)
            return np.stack([s(x % 1 if periodic else x) for s in splines],
                            axis=1)

        return UnivariateFunction(f, periodic=periodic)

    else:
        c0_indices = pc.c0_indices(thresh=c0_thresh)

        if len(c0_indices) == 0:
            return InterpolatedUnivariateSpline(pc,
                                                c0_thresh=None,
                                                k=k,
                                                **scipyargs)

        shift_amount = 0

        if pc.periodic:  # roll pointcloud to first C0-index, then it can be treated as an unperiodic Pointcloud
            shift_amount = pc.verts[c0_indices[0]]
            pc = pc.roll(-c0_indices[0])
            c0_indices = (
                c0_indices -
                c0_indices[0])[1:]  # first index is gonna be added anyways

        c0_indices = [0] + c0_indices.tolist() + [
            pc._points.shape[0] - 1
        ]  # pc._points because in the periodic case we repeat the last

        funcs = [
            InterpolatedUnivariateSpline(pc._points[i:j + 1],
                                         c0_thresh=None,
                                         k=k,
                                         **scipyargs)
            for i, j in zip(c0_indices[:-1], c0_indices[1:])
        ]

        return UnivariateFunction.stack_multiple(
            funcs, pc._verts[c0_indices[1:-1]],
            periodic=pc.periodic).shift(shift_amount)
Ejemplo n.º 10
0
 def toPointCloud(self, x, verts=None):
     if self._periodic:
         if x[-1] != x[0]:
             log.warning( 'Warning, a periodic function is evaluated over' + \
                 ' a non-periodic set of points, this is usually a bug' )
     return PointCloud(self(x), verts=verts)
Ejemplo n.º 11
0
def main(
    nelems: 'number of elements' = 20,
    epsilon: 'epsilon, 0 for automatic (based on nelems)' = 0,
    timestep: 'time step' = .01,
    maxtime: 'end time' = 1.,
    theta: 'contact angle (degrees)' = 90,
    init: 'initial condition (random/bubbles)' = 'random',
    withplots: 'create plots' = True,
):

    mineps = 1. / nelems
    if not epsilon:
        log.info('setting epsilon=%f' % mineps)
        epsilon = mineps
    elif epsilon < mineps:
        log.warning('epsilon under crititical threshold: %f < %f' %
                    (epsilon, mineps))

    ewall = .5 * numpy.cos(theta * numpy.pi / 180)

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

    # prepare bases
    cbasis, mubasis = function.chain(
        [domain.basis('spline', degree=2),
         domain.basis('spline', degree=2)])

    # define mixing energy and splitting: F' = f_p - f_n
    F = lambda c_: (.5 / epsilon**2) * (c_**2 - 1)**2
    f_p = lambda c_: (1. / epsilon**2) * 4 * c_
    f_n = lambda c_: (1. / epsilon**2) * (6 * c_ - 2 * c_**3)

    # prepare matrix
    A = function.outer( cbasis ) \
      + (timestep*epsilon**2) * function.outer( cbasis.grad(geom), mubasis.grad(geom) ).sum(-1) \
      + function.outer( mubasis, mubasis - f_p(cbasis) ) \
      - function.outer( mubasis.grad(geom), cbasis.grad(geom) ).sum(-1)
    matrix = domain.integrate(A, geometry=geom, ischeme='gauss4')

    # prepare wall energy right hand side
    rhs0 = domain.boundary.integrate(mubasis * ewall,
                                     geometry=geom,
                                     ischeme='gauss4')

    # construct initial condition
    if init == 'random':
        numpy.random.seed(0)
        c = cbasis.dot(numpy.random.normal(0, .5, cbasis.shape))
    elif init == 'bubbles':
        R1 = .25
        R2 = numpy.sqrt(.5) * R1  # area2 = .5 * area1
        c = 1 + function.tanh( (R1-function.norm2(geom-(.5+R2/numpy.sqrt(2)+.8*epsilon)))/epsilon ) \
              + function.tanh( (R2-function.norm2(geom-(.5-R1/numpy.sqrt(2)-.8*epsilon)))/epsilon )
    else:
        raise Exception('unknown init %r' % init)

    # prepare plotting
    nsteps = numeric.round(maxtime / timestep)
    makeplots = MakePlots(domain, geom, nsteps, video=withplots
                          == 'video') if withplots else lambda *args: None

    # start time stepping
    for istep in log.range('timestep', nsteps):

        Emix = F(c)
        Eiface = .5 * (c.grad(geom)**2).sum(-1)
        Ewall = (abs(ewall) + ewall * c)

        b = cbasis * c - mubasis * f_n(c)
        rhs, total, energy_mix, energy_iface = domain.integrate(
            [b, c, Emix, Eiface], geometry=geom, ischeme='gauss4')
        energy_wall = domain.boundary.integrate(Ewall,
                                                geometry=geom,
                                                ischeme='gauss4')
        log.user('concentration {}, energy {}'.format(
            total, energy_mix + energy_iface + energy_wall))

        lhs = matrix.solve(rhs0 + rhs, tol=1e-12, restart=999)
        c = cbasis.dot(lhs)
        mu = mubasis.dot(lhs)

        makeplots(c, mu, energy_mix, energy_iface, energy_wall)

    return lhs, energy_mix, energy_iface, energy_wall
Ejemplo n.º 12
0
def main(
    nelems: 'number of elements' = 20,
    epsilon: 'epsilon, 0 for automatic (based on nelems)' = 0,
    timestep: 'time step' = .01,
    maxtime: 'end time' = 1.,
    theta: 'contact angle (degrees)' = 90,
    init: 'initial condition (random/bubbles)' = 'random',
    withplots: 'create plots' = True,
  ):

  mineps = 1./nelems
  if not epsilon:
    log.info('setting epsilon={}'.format(mineps))
    epsilon = mineps
  elif epsilon < mineps:
    log.warning('epsilon under crititical threshold: {} < {}'.format(epsilon, mineps))

  # create namespace
  ns = function.Namespace()
  ns.epsilon = epsilon
  ns.ewall = .5 * numpy.cos( theta * numpy.pi / 180 )

  # construct mesh
  xnodes = ynodes = numpy.linspace(0,1,nelems+1)
  domain, ns.x = mesh.rectilinear( [ xnodes, ynodes ] )

  # prepare bases
  ns.cbasis, ns.mubasis = function.chain([
    domain.basis('spline', degree=2),
    domain.basis('spline', degree=2)
  ])

  # polulate namespace
  ns.c = 'cbasis_n ?lhs_n'
  ns.c0 = 'cbasis_n ?lhs0_n'
  ns.mu = 'mubasis_n ?lhs_n'
  ns.f = '(6 c0 - 2 c0^3 - 4 c) / epsilon^2'

  # construct initial condition
  if init == 'random':
    numpy.random.seed(0)
    lhs0 = numpy.random.normal(0, .5, ns.cbasis.shape)
  elif init == 'bubbles':
    R1 = .25
    R2 = numpy.sqrt(.5) * R1 # area2 = .5 * area1
    ns.cbubble1 = function.tanh((R1-function.norm2(ns.x-(.5+R2/numpy.sqrt(2)+.8*ns.epsilon)))/ns.epsilon)
    ns.cbubble2 = function.tanh((R2-function.norm2(ns.x-(.5-R1/numpy.sqrt(2)-.8*ns.epsilon)))/ns.epsilon)
    sqr = domain.integral('(c - cbubble1 - cbubble2 - 1)^2 + mu^2' @ ns, geometry=ns.x, degree=4)
    lhs0 = solver.optimize('lhs', sqr)
  else:
    raise Exception( 'unknown init %r' % init )

  # construct residual
  res = domain.integral('epsilon^2 cbasis_n,k mu_,k + mubasis_n (mu + f) - mubasis_n,k c_,k' @ ns, geometry=ns.x, degree=4)
  res -= domain.boundary.integral('mubasis_n ewall' @ ns, geometry=ns.x, degree=4)
  inertia = domain.integral('cbasis_n c' @ ns, geometry=ns.x, degree=4)

  # solve time dependent problem
  nsteps = numeric.round(maxtime/timestep)
  makeplots = MakePlots(domain, nsteps, ns) if withplots else lambda *args: None
  for istep, lhs in log.enumerate('timestep', solver.impliciteuler('lhs', target0='lhs0', residual=res, inertia=inertia, timestep=timestep, lhs0=lhs0)):
    makeplots(lhs)
    if istep == nsteps:
      break

  return lhs0, lhs