def av_integrate(avfun, avmap, N):
    """Approximate the integral of a function of active variables.

    Parameters
    ----------
    avfun : function
        a function of the active variables
    avmap : ActiveVariableMap
        a domains.ActiveVariableMap
    N : int
        the number of points in the quadrature rule

    Returns
    -------
    mu : float
        an estimate of the integral

    Notes
    -----
    This function is usually used when one has already constructed a response
    surface on the active variables and wants to estimate its integral.
    """
    if not isinstance(avmap, ActiveVariableMap):
        raise TypeError('avmap should be an ActiveVariableMap.')

    if not isinstance(N, Integral):
        raise TypeError('N should be an integer.')

    Yp, Yw = av_quadrature_rule(avmap, N)
    if isinstance(avfun, ActiveSubspaceResponseSurface):
        avf = avfun.predict_av(Yp)[0]
    else:
        avf = SimulationRunner(avfun).run(Yp)
    mu = np.dot(Yw.T, avf)[0, 0]
    return mu
예제 #2
0
    def build_from_interface(self, m, fun, dfun=None):
        self.m = m

        # number of gradient samples
        M = 6 * (m + 1) * np.log(m)

        # sample points
        if self.bflag:
            X = np.random.uniform(-1.0, 1.0, size=(M, m))
        else:
            X = np.random.normal(size=(M, m))
        self.X = X

        # sample the simulation's outputs
        sr = SimulationRunner(fun)
        self.f_runner = sr

        # sample the simulation's gradients
        if dfun == None:
            df = finite_difference_gradients(X, sr)
        else:
            sgr = SimulationGradientRunner(dfun)
            df = sgr.run(X)
            self.df_runner = sgr
        self.df = df

        # compute the active subspaces
        self.subspaces = build_subspaces(df)
예제 #3
0
def finite_difference_gradients(X, fun, h=1e-6):
    """Compute finite difference gradients with a given interface.

    Parameters
    ----------
    X : ndarray 
        M-by-m matrix that contains the points to estimate the gradients with 
        finite differences
    fun : function
        function that returns the simulation's quantity of interest given inputs
    h : float, optional 
        the finite difference step size (default 1e-6)

    Returns
    -------
    df : ndarray 
        M-by-m matrix that contains estimated partial derivatives approximated 
        by finite differences
    """
    X, M, m = process_inputs(X)

    # points to run simulations including the perturbed inputs
    XX = np.kron(np.ones((m+1, 1)),X) + \
        h*np.kron(np.vstack((np.zeros((1, m)), np.eye(m))), np.ones((M, 1)))

    # run the simulation
    if isinstance(fun, SimulationRunner):
        F = fun.run(XX)
    else:
        F = SimulationRunner(fun).run(XX)

    df = (F[M:].reshape((m, M)).transpose() - F[:M]) / h
    return df.reshape((M,m))
 def train_with_interface(self, fun, N, NMC=10):
     Y, X, ind = as_design(self.avmap, N, NMC=NMC)
     
     if isinstance(self.avmap.domain, BoundedActiveVariableDomain):
         X = np.vstack((X, self.avmap.domain.vertX))
         Y = np.vstack((Y, self.avmap.domain.vertY))
         il = np.amax(ind) + 1
         iu = np.amax(ind) + self.avmap.domain.vertX.shape[0] + 1
         iind = np.arange(il, iu)
         ind = np.vstack(( ind, iind.reshape((iind.size,1)) ))
         
     # run simulation interface at all design points
     sr = SimulationRunner(fun)
     f = sr.run(X)
     
     Ef, Vf = conditional_expectations(f, ind)
     self.train(Y, Ef, v=Vf)
예제 #5
0
    def build_from_interface(self, m, fun, dfun=None, avdim=None):
        self.m = m

        # number of gradient samples
        M = int(np.floor(6*(m+1)*np.log(m)))

        # sample points for gradients
        if self.bndflag:
            X = np.random.uniform(-1.0, 1.0, size=(M, m))
        else:
            X = np.random.normal(size=(M, m))
        funr = SimulationRunner(fun) 
        f = funr.run(X)
        self.X, self.f, self.funr = X, f, funr
        
        # sample the simulation's gradients
        if dfun == None:
            df = finite_difference_gradients(X, fun)
        else:
            dfunr = SimulationGradientRunner(dfun)
            df = dfunr.run(X)
            self.dfunr = dfunr

        # compute the active subspace
        ss = Subspaces()
        ss.compute(df)
        if avdim is not None:
            ss.partition(avdim)
        self.n = ss.W1.shape[1]
        print 'Dimension of subspace is {:d}'.format(self.n)
        
        # set up the active variable domain and map
        if self.bndflag:
            avdom = BoundedActiveVariableDomain(ss)
            avmap = BoundedActiveVariableMap(avdom)
        else:
            avdom = UnboundedActiveVariableDomain(ss)
            avmap = UnboundedActiveVariableMap(avdom)
            
        # build the response surface
        avrs = ActiveSubspaceResponseSurface(avmap)
        avrs.train_with_interface(fun, int(np.power(5,self.n)))
        self.av_respsurf = avrs
    def train_with_interface(self, fun, N, NMC=10):
        """Train the response surface with input/output pairs.

        Parameters
        ----------
        fun : function 
            a function that returns the simulation quantity of interest given a 
            point in the input space as an 1-by-m ndarray
        N : int 
            the number of points used in the design-of-experiments for 
            constructing the response surface
        NMC : int, optional 
            the number of points used to estimate the conditional expectation 
            and conditional variance of the function given a value of the active
            variables

        Notes
        -----
        The training methods exploit the eigenvalues from the active subspace
        analysis to determine length scales for each variable when tuning
        the parameters of the radial bases.

        The method sets attributes of the object for further use.

        The method uses the response_surfaces.av_design function to get the
        design for the appropriate `avmap`.
        """
        Y, X, ind = av_design(self.avmap, N, NMC=NMC)

        if isinstance(self.avmap.domain, BoundedActiveVariableDomain):
            X = np.vstack((X, self.avmap.domain.vertX))
            Y = np.vstack((Y, self.avmap.domain.vertY))
            il = np.amax(ind) + 1
            iu = np.amax(ind) + self.avmap.domain.vertX.shape[0] + 1
            iind = np.arange(il, iu)
            ind = np.vstack((ind, iind.reshape((iind.size, 1))))

        # run simulation interface at all design points
        if isinstance(fun, SimulationRunner):
            f = fun.run(X)
        else:
            f = SimulationRunner(fun).run(X)

        Ef, Vf = conditional_expectations(f, ind)
        self._train(Y, Ef, v=Vf)
def as_integrate(fun, domain, subspace, N, NMC=10):
    """
    Description of as_integrate

    Arguments:
        fun:
        domain:
        subspace:
        N:
        NMC: (default=10)
    Outputs:
        mu:
        lb:
        ub:
    """
    w, x, ind = as_quadrature_rule(domain, subspace, N, NMC)
    f = SimulationRunner(fun).run(x)
    return compute_integral(w, f, ind)
예제 #8
0
    def build_from_interface(self, fun, dfun=None, avdim=None):
        """
        Build the active subspace-enabled model with interfaces to the
        simulation.

        :param function fun: A function that interfaces with the simulation.
            It should take an ndarray of shape 1-by-m (e.g., a row of `X`), and
            it should return a scalar. That scalar is the quantity of interest from the simulation.
        :param function dfun: A function that interfaces with the simulation.
            It should take an ndarray of shape 1-by-m (e.g., a row of `X`), and it
            should return the gradient of the quantity of interest as an ndarray of shape 1-by-m.
        :param int avdim: The dimension of the active subspace. If `avdim` is not
            present, it is chosen after computing the active subspaces using
            the given interfaces.

        **Notes**

        This method follows these steps:

        #. Draw random points according to the weight function on the space\
        of simulation inputs.
        #. Compute the quantity of interest and its gradient at the sampled\
        inputs. If `dfun` is None, use finite differences.
        #. Use the collection of gradients to estimate the eigenvectors and\
        eigenvalues that determine and define the active subspace.
        #. Train a response surface using the interface, which uses a careful\
        design of experiments on the space of active variables. This design\
        uses about 5 points per dimension of the active subspace.

        """
        if not hasattr(fun, '__call__'):
            raise TypeError('fun should be a callable function.')

        if dfun is not None:
            if not hasattr(dfun, '__call__'):
                raise TypeError('dfun should be a callable function.')

        if avdim is not None:
            if not isinstance(avdim, int):
                raise TypeError('avdim should be an integer')

        m = self.m

        # number of gradient samples
        M = int(np.floor(6*(m+1)*np.log(m)))

        # sample points for gradients
        if self.bounded_inputs:
            X = np.random.uniform(-1.0, 1.0, size=(M, m))
        else:
            X = np.random.normal(size=(M, m))

        fun = SimulationRunner(fun)
        f = fun.run(X)
        self.X, self.f, self.fun = X, f, fun

        # sample the simulation's gradients
        if dfun == None:
            df = finite_difference_gradients(X, fun)
        else:
            dfun = SimulationGradientRunner(dfun)
            df = dfun.run(X)
            self.dfun = dfun

        # compute the active subspace
        ss = Subspaces()
        ss.compute(df)

        if avdim is not None:
            ss.partition(avdim)
        self.n = ss.W1.shape[1]
        print 'The dimension of the active subspace is {:d}.'.format(self.n)

        # set up the active variable domain and map
        if self.bounded_inputs:
            avdom = BoundedActiveVariableDomain(ss)
            avmap = BoundedActiveVariableMap(avdom)
        else:
            avdom = UnboundedActiveVariableDomain(ss)
            avmap = UnboundedActiveVariableMap(avdom)

        # build the response surface
        asrs = ActiveSubspaceResponseSurface(avmap)
        asrs.train_with_interface(fun, int(np.power(5,self.n)))

        # compute testing error as an R-squared
        self.Rsqr = 1.0 - ( np.linalg.norm(asrs.predict(X)[0] - f)**2 \
                            / np.var(f) )

        self.as_respsurf = asrs
def integrate(fun, avmap, N, NMC=10):
    """Approximate the integral of a function of m variables.

    Parameters
    ----------
    fun : function
        an interface to the simulation that returns the quantity of interest
        given inputs as an 1-by-m ndarray
    avmap : ActiveVariableMap
        a domains.ActiveVariableMap
    N : int
        the number of points in the quadrature rule
    NMC : int, optional
        the number of points in the Monte Carlo estimates of the conditional
        expectation and conditional variance (default 10)

    Returns
    -------
    mu : float
        an estimate of the integral of the function computed against the weight
        function on the simulation inputs
    lb : float
        a central-limit-theorem 95% lower confidence from the Monte Carlo part
        of the integration
    ub : float
        a central-limit-theorem 95% upper confidence from the Monte Carlo part
        of the integration

    See Also
    --------
    integrals.quadrature_rule

    Notes
    -----
    The CLT-based bounds `lb` and `ub` are likely poor estimators of the error.
    They only account for the variance from the Monte Carlo portion. They do
    not include any error from the integration rule on the active variables.
    """
    if not isinstance(avmap, ActiveVariableMap):
        raise TypeError('avmap should be an ActiveVariableMap.')

    if not isinstance(N, Integral):
        raise TypeError('N should be an integer')

    # get the quadrature rule
    Xp, Xw, ind = quadrature_rule(avmap, N, NMC=NMC)

    # compute the simulation output at each quadrature node
    if isinstance(fun, SimulationRunner):
        f = fun.run(Xp)
    else:
        f = SimulationRunner(fun).run(Xp)

    # estimate conditional expectations and variances
    Ef, Vf = conditional_expectations(f, ind)

    # get weights for the conditional expectations
    w = conditional_expectations(Xw * NMC, ind)[0]

    # estimate the average
    mu = np.dot(Ef.T, w)

    # estimate the variance due to Monte Carlo
    sig2 = np.dot(Vf.T, w * w) / NMC

    # compute 95% confidence bounds from the Monte Carlo
    lb, ub = mu - 1.96 * np.sqrt(sig2), mu + 1.96 * np.sqrt(sig2)
    return mu[0, 0], lb[0, 0], ub[0, 0]