iam-git / WellMet (public) (License: MIT) (since 2021-08-31) (hash sha1)
WellMet is pure Python framework for spatial structural reliability analysis. Or, more specifically, for "failure probability estimation and detection of failure surfaces by adaptive sequential decomposition of the design domain".

/convex_hull.py (5106dbfd64b5403599d51598c03d9fad60956286) (26963 bytes) (mode 100644) (type blob)

#!/usr/bin/env python
# coding: utf-8

import numpy as np
import mpmath
from scipy import stats
from scipy import spatial # for QHull

mpmath.mp.dps = 325 # to cover anything that double-presigion float can handle

# maximum radius, where norm.pdf() wasn't zero
# -38.575500173381374935388521407730877399444580.....
# don't ask me what the magic python use to distinguish
# digits after double precision
max_R_ever = 38.5755

#č jako sůl potřebujem statefull třidu,
#č která by schovávala vnitřní implementaciju,
#č která by nabízela jednotné rozhraní pro ostatní kód WellMet
#č a byla by přiměřeně kompatibilní s ConvexHull ze scipy.spatial
#č Ze scipy bych viděl podporu atributů .points, .npoints, .equations
#č Určitě bych chtěl funkce .update() a .is_outside(sample)
#č Atribute .space by ukazoval, v jakém prostoru konvexní obál je vytvořen
# GBall, BrickHull, DirectHull, CompleteHull and QHull has been implemented


def calculate_brick_complement_probability(mins, maxs):
    "For Gaussian space"
    #č na začátku nastavíme obyčejnou 1
    #č v cyklu bude nasobit delkami 
    #č podle jednotlivých souřadnic
    #č a typ se změní na mpf (hodně přesná aritmetika)
    volume = 1
    
    # ncdf() do not support matrix input
    for xmin, xmax in zip(mins, maxs):
        volume *= mpmath.ncdf(xmax) - mpmath.ncdf(xmin)
    
    # print test
#    if len(mins) == 2:
#        min_x, min_y = mins
#        max_x, max_y = maxs
#        test_volume = 0
#        test_volume += np.sum(stats.norm.cdf(mins)) 
#        test_volume += np.sum(stats.norm.sf(maxs))
#        test_volume -= stats.norm.cdf(min_x) * stats.norm.cdf(min_y)
#        test_volume -= stats.norm.cdf(min_x) * stats.norm.sf(max_y)
#        test_volume -= stats.norm.sf(max_x) * stats.norm.cdf(min_y)
#        test_volume -= stats.norm.sf(max_x) * stats.norm.sf(max_y)
#        print("test_volume:", test_volume)  
        
    return float(1-volume)
    
# inspired by 
# https://math.stackexchange.com/questions/192610/how-to-build-a-orthogonal-basis-from-a-vector/979013
# most interesting from there:
# """Pick a vector. WLOG, you chose $(x_1,x_2,x_3,x_4)$. 
# Now write it as a quaternion: $$x_1+ix_2+jx_3+kx_4$$ 
# Then, since multiplication by $i,j,k$ 
# rotates this vector $90^0$ across the various axes of our 4D space,
# the following three vectors make your initial choice of vector into an orthonormal basis: 
# $$i(x_1+ix_2+jx_3+kx_4)=ix_1-x_2+kx_3-jx_4\mapsto (-x_2,x_1,-x_4,x_3)$$ 
# $$j(x_1+ix_2+jx_3+kx_4)=jx_1-kx_2-x_3+ix_4\mapsto (-x_3,x_4,x_1,-x_2)$$
# $$k(x_1+ix_2+jx_3+kx_4)=kx_1+jx_2-ix_3-x_4\mapsto (-x_4,-x_3,x_2,x_1)$$"""
def get_orth_basis(vector):
    dim = len(vector)
    if dim == 2:
        x, y = vector
        return np.array([[x, y],[-y, x]])
    elif dim == 4:
        x_1, x_2, x_3, x_4 = vector
        return np.array([
            [x_1, x_2, x_3, x_4],
            [-x_2, x_1, -x_4, x_3],
            [-x_3, x_4, x_1, -x_2],
            [-x_4, -x_3, x_2, x_1]])
    else:
        random_basis = np.random.random((dim, dim))
        #č QR rozklad jede po sloupcich
        random_basis[:,0] = vector
        #č co jsem viděl, numpy matici Q normalizuje
        #č a první sloupec zůstavá (skoro) tím samým, co byl před tím
        Q, __R = np.linalg.qr(random_basis)
        # yep. Numpy sometimes inverts sign of the first Q matrix vector
        Q[:,0] = vector
        return Q.T
    
        
        
#ё рельса
#č nepodařílo se mi nějak rozumně zobecnit pro libovolný prostor
#č takže Gauss
#č (každá třida implementuje zvlášť)

def shot(hull, ns, use_MC=False):
    if hull.space == 'G':
        to_fire = np.nanargmax(hull.b)
        r = -hull.b[to_fire]
        a = hull.A[to_fire]
        
        if use_MC:
            fire_from = stats.norm.sf(r)
            t = np.linspace(fire_from, 0, ns, endpoint=False)
            t = stats.norm.isf(t)
        else:
            if r < max_R_ever:
                R = max_R_ever
            else:
                R = r + 10
            t = np.linspace(r, R, ns, endpoint=True)
            
        fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
        return hull.sample.f_model.new_sample(fire_G, space='G')


def fire(hull, ns, use_MC=False):
    if hull.space == 'G':
        to_fire = np.nanargmax(hull.b)
        r = -hull.b[to_fire]
        a = hull.A[to_fire]
        
        orth_nodes_T = np.random.randn(len(a), ns) # len(a) == ndim
        orth_basis = get_orth_basis(a)
        
        if use_MC:
            fire_from = stats.norm.sf(r)
            t = np.linspace(fire_from, 0, ns, endpoint=False)
            orth_nodes_T[0] = stats.norm.isf(t)
        else:
            if r < max_R_ever:
                R = max_R_ever
            else:
                R = r + 10
            orth_nodes_T[0] = np.linspace(r, R, ns, endpoint=True)
            
        fire_G = (orth_basis.T @ orth_nodes_T).T
        return hull.sample.f_model.new_sample(fire_G, space='G')

#
# Common orth functions for DirectHull, CompleteHull and QHull
#

# (even private function can be shared :D )
def _orth_helper(hull):
    # supposed hull.space == 'G'
    hull._update()
    to_fire = np.nanargmax(hull.b)
    a = hull.A[to_fire]
    
    orth_basis = get_orth_basis(a)
    
    #č musí tam bejt G coordinates
    A = orth_basis @ hull.points.T
    bp = np.nanmax(A, axis=1)
    bm = np.nanmin(A, axis=1)
    
    return orth_basis, bp, bm
    
    
def get_orth_outside(hull):
    if hull.space == 'G' :
        # hull._update() will perform _orth_helper()
        A, bp, bm = hull._orth_helper()
        return calculate_brick_complement_probability(bm, bp)
    else:
        #č když prostor není G, můžeme sice něco odvodit od s-ball 
        #č ale nechť se s tím pará GHull, vrátíme nulu.
        return 0

def get_orth_equations(hull):
    # hull._update() will perform _orth_helper()
    direct_plan, bp, bm = hull._orth_helper()

    A = np.vstack((-direct_plan, direct_plan))
    b = np.concatenate((bm, -bp))
    return np.hstack((A,b[:,None]))     


#
# Common 2FORM functions for DirectHull, CompleteHull and QHull
#

# (even private function can be shared :D )
def _2FORM_helper(hull):
    # supposed hull.space == 'G'
    hull._update()
    to_fire = np.nanargmax(hull.b)
    a = hull.A[to_fire]
    
    #č musí tam bejt G coordinates
    # we'll get 1D npoints-sized negative b array
    x = a @ hull.points.T
    # here x == _b array has another meaning than hull.b
    # where in hull.b we search for least distanced hyperplane
    # here, in x, we search for a valid b_backward hyperplane offset
    # It can be seen as transformation of the ED to 1D a-projection.
    # There evidently exist min(x) and max(x)
    # such that every single x value lies inside (between)
    # max(x) (forward) value should be corresponding to r, to -hull.b
    # So, we need to find min(x) (backward) value
    # x_forward == max(x) == -b == r == max(a @ points.T) == -max(hull.b) == min(-hull.b)
    # x_backward == min(x) == max(-a @ points.T)
    x_forward = np.max(x)
    x_backward = np.min(x)
    #print(-np.max(hull.b), x_forward) 
    return x_backward, x_forward, a


def get_2FORM_outside(hull):
    if hull.space == 'G':
        # hull._update() will perform _FORM_helper()
        x_backward, x_forward, __a = hull._2FORM_helper()
        return stats.norm.cdf(x_backward) + stats.norm.sf(x_forward)
    else:
        #č když prostor není G, můžeme sice něco odvodit od s-ball 
        #č ale nechť se s tím pará GHull, vrátíme nulu.
        return 0
            
def get_2FORM_equations(hull):
    # hull._update() will perform _2FORM_helper()
    x_backward, x_forward, a = hull._2FORM_helper()
    A = np.vstack((-a, a))
    b = np.array([[x_backward],[-x_forward]])
    return np.hstack((A,b))  
    

#č jistě musíme mít nějaký zbytečný kus kódu
#č třida jen pro formu, jen tak na hračku
#č když tečíčky jsou dále jak nejvzdálenější vzorek (bod),
#č tak řekneme, že jsou totálně mimo!
class GBall:
    def __init__(hull, sample):
        hull.sample = sample
        hull.space = 'G' #оӵ малы? Кинлы со кулэ?
    
    @property
    def points(hull):
        return hull.sample.G
        
    @property
    def npoints(hull):
        return len(hull.sample)
        
    @property
    def nsimplex(hull):
        return 1
        
    #def update(hull): pass
    
    def is_inside(hull, nodes):
        R2_hull = np.nanmax(np.sum(np.square(hull.sample.G), axis=1))
        
        R2_nodes = np.sum(np.square(nodes.G), axis=1)
        return R2_nodes < R2_hull 

    def is_outside(hull, nodes): 
        return ~hull.is_inside(nodes)
        
    def get_R(hull):
        return np.sqrt(np.nanmax(np.sum(np.square(hull.sample.G), axis=1)))
        
    def get_r(hull): 
        # zero at best, 
        # on the safe side, between -R and 0
        #r2_hull = np.nanmin(np.sum(np.square(hull.sample.G), axis=1))
        #return -np.sqrt(r2_hull)
        
        #č když na get_r nahlížím z hlediska toho,
        #č že r vymezuje mezikruží, kde nejsme jistí
        #č inside-outside,
        #č tak zde get_r vždy musí vracet prostě R
        return hull.get_R()
        
    def get_orth_outside(hull):
        x = np.full(hull.sample.nvar, hull.get_r())
        return calculate_brick_complement_probability(-x, x)
        
     
    def get_2FORM_outside(hull):
        return 2*stats.norm.sf(hull.get_r())
        
     
    def fire(hull, *args, **kwargs):
        pass
    
    def shot(hull, *args, **kwargs):
        pass
        


#č vyhejbám slovu Box, ať nevnáším ještě víc zmatku v označení
class BrickHull: #č nebo BoundingBrick
    def __init__(hull, sample, space='G'):
        hull.sample = sample
        hull.space = space
        
        
        hull._npoints = 0
        #č miny a maxy obsahují minima a maxima 
        #č podel jednotlivých souřadnic (proměnných)
        #č sincerely yours, Captain Obvious.
        hull.mins = np.full(hull.sample.nvar, np.inf)
        hull.maxs = np.full(hull.sample.nvar, -np.inf)
    
    def _update(hull):
        if hull._npoints < hull.npoints:
            hull.mins = np.nanmin(hull.points, axis=0)
            hull.maxs = np.nanmax(hull.points, axis=0)
            hull._npoints = hull.npoints
        
        
    def is_inside(hull, nodes):
        hull._update()
        x = getattr(nodes, hull.space)
        more = np.all(x < hull.maxs, axis=1)
        less = np.all(x > hull.mins, axis=1)
        
        return np.all((more, less), axis=0) #np.all(np.array((more, less)), axis=0)
        
    def is_outside(hull, nodes): 
        return ~hull.is_inside(nodes)
        
    def get_design_points(hull):
        hull._update()
        #sample_model = -hull.A * hull.b.reshape(-1,1)
        sample_model = np.vstack((np.diag(hull.maxs), np.diag(hull.mins)))
        return hull.sample.f_model.new_sample(sample_model, space=hull.space)
        
    # shortcut for Ghull
    # valid only if space==G
    def get_r(hull):
        if hull.space=='G':
            hull._update()
            return np.min((-hull.mins, hull.maxs))
        else:
            return 0
        
    @property
    def points(hull):
        return getattr(hull.sample, hull.space)
        
    @property
    def npoints(hull):
        return len(hull.sample)
        
    @property
    def nsimplex(hull):
        return hull.sample.nvar*2
        
    @property
    def A(hull):
        hull._update()
        #č žádná optimizace, ale kdo tu funkci bude spouštět?
        diag = np.diag(np.ones(hull.sample.nvar))
        return np.vstack((diag, -diag))
        
    @property
    def b(hull):
        hull._update()
        return np.concatenate((hull.maxs, -hull.mins))
        
    @property
    def equations(hull):
        hull._update()
        #č žádná optimizace, ale kdo tu funkci bude spouštět?
        diag = np.diag(np.ones(hull.sample.nvar))
        
        A = np.vstack((diag, -diag))
        b = np.concatenate((-hull.maxs, hull.mins))
        return np.hstack((A,b[:,None]))
        
        
    def get_orth_outside(hull):
        if hull.space == 'G':
            hull._update()
            return calculate_brick_complement_probability(hull.mins, hull.maxs)
        else:
            #č když prostor není G, můžeme sice něco odvodit od s-ball 
            #č ale nechť se s tím pará GHull, vrátíme nulu.
            return 0
            
    def get_orth_equations(hull):
        return hull.equations
        
        
    def _2FORM_helper(hull):
        hull._update()
        #č je to úplně zbýtečně trapit se nejakejma argama
        #č poďme na to globálně!
        #č (doufám, že nvar*2 výpočtů cdf() nikomu neublíží... )
        #min_arg = np.nanargmin(hull.mins)
        #max_arg = np.nanargmax(hull.maxs)
        cdfs = stats.norm.cdf(hull.mins)
        sfs = stats.norm.sf(hull.maxs)
        #č _update() počítá nanmin a nanmax
        #č takže nemusíme se bat nějakého hnusu v číslech
        pfs = cdfs + sfs
        i = np.argmax(pfs) # number of variable
        # we'll return 2FORM estimation 
        # and corresponding number of variable
        return pfs[i], i
        
        
    def get_2FORM_outside(hull):
        if hull.space == 'G':
            # hull._update() will perform _FORM_helper()
            p_out, __i = hull._2FORM_helper()
            return p_out
        else:
            #č když prostor není G, můžeme sice něco odvodit od s-ball 
            #č ale nechť se s tím pará GHull, vrátíme nulu.
            return 0
            
        
    def get_2FORM_equations(hull):
        __p_out, i = hull._2FORM_helper()
        equations = np.zeros((2,hull.sample.nvar+1))
        equations[0, i] = 1
        equations[0, -1] = -hull.maxs[i]
        equations[1, i] = -1
        equations[1, -1] = hull.mins[i]
        return equations
        
    def shot(hull, *args, **kwargs):
        pass
        
        # add use_MC to the function signature
        # but remain unimplemented for now
        #č ten fire je hrozný. Výhodit U prostor 
        #č (jsou tam nepřesné cdf transformace)
        #č po úpravě přídat hull._update() na začátku.
    def fire(hull, ns, use_MC=False): # boom
        sample_U = hull.sample.U
        
        U_mins = np.nanmin(sample_U, axis=0)
        U_maxs = np.nanmax(sample_U, axis=0)
        
        #č min nebo max?
        if np.nanmax(U_mins) > (1-np.nanmin(U_maxs)):
            #č expandujeme se dolu
            var = np.nanargmax(U_mins)
            value = U_mins[var]
            t = np.linspace(value, 0, ns, endpoint=False)
        else:
            #č expandujeme se nahoru
            var = np.nanargmin(U_maxs)
            value = U_maxs[var]
            t = np.linspace(value, 1, ns, endpoint=False)
        
        boom = np.random.random((ns, hull.sample.nvar))
        boom[:, var] = t
        
        return hull.sample.f_model.new_sample(boom, space='U')
    

class DirectHull:
    # take some global functions
    fire = fire
    shot = shot
    _orth_helper = _orth_helper
    get_orth_outside = get_orth_outside
    get_orth_equations = get_orth_equations
    
    _2FORM_helper = _2FORM_helper
    get_2FORM_outside = get_2FORM_outside
    get_2FORM_equations = get_2FORM_equations
    
    
    def __init__(hull, sample, direct_plan, space='G'):
        hull.sample = sample
        hull.direct_plan = direct_plan
        hull.space = space
        
        hull.regen()
        
    def regen(hull):
        hull._npoints = 0
        
        hull.bp = np.full(len(hull.direct_plan), -np.inf)
        #hull.update()
    
    #č nejsem jist, jestli ten update vůbec dělat.
    #č lze navrhnout třidu aj tak, že sama bude hlídat změny.
    #č Jenomže, co kdybychom ten automatický update nechtěli?
    def _update(hull):
        if hull._npoints < hull.npoints:
            #hull.points = getattr(hull.sample, hull.space)
            new_points = hull.points[hull._npoints:]
            
            A = hull.direct_plan @ new_points.T
            new_bp = np.nanmax(A, axis=1)
            hull.bp = np.nanmax((new_bp, hull.bp), axis=0)
            
            hull._npoints = hull.npoints

    @property
    def points(hull):
        return getattr(hull.sample, hull.space)
        
    @property
    def npoints(hull):
        return len(hull.sample)
        
    @property
    def nsimplex(hull):
        return len(hull.direct_plan)
        
    @property
    def A(hull):
        return hull.direct_plan
        
    @property
    def b(hull):
        hull._update()
        return -hull.bp
        
    @property
    def equations(hull):
        return np.hstack((hull.A, hull.b[:,None]))

    def is_inside(hull, nodes):
        hull._update()
        x = getattr(nodes, hull.space)
        
        #E [normal, offset] forming the hyperplane equation of the facet (see Qhull documentation for more)
        A = hull.direct_plan
        bp = np.atleast_2d(hull.bp).T
        
        # N=ns, E - number of hyperplane equations
        ExN = A @ x.T
        higher = np.all(ExN < bp, axis=0)
        return higher
    
    def is_outside(hull, nodes): 
        return ~hull.is_inside(nodes)
        
    def get_design_points(hull):
        hull._update()
        sample_model = -hull.A * hull.b.reshape(-1,1)
        return hull.sample.f_model.new_sample(sample_model, space=hull.space)
        
    # shortcut for Ghull
    # valid only if space==G
    def get_r(hull):
        if hull.space=='G':
            hull._update()
            return np.min(hull.bp)
        else:
            return 0
    
    


class CompleteHull:
    # take some global functions
    fire = fire
    shot = shot
    _orth_helper = _orth_helper
    get_orth_outside = get_orth_outside
    get_orth_equations = get_orth_equations
    
    _2FORM_helper = _2FORM_helper
    get_2FORM_outside = get_2FORM_outside
    get_2FORM_equations = get_2FORM_equations
    
    def __init__(hull, sample, direct_plan, space='G'):
        hull.sample = sample
        hull.direct_plan = direct_plan
        hull.space = space
        
        hull.regen()
        
    def regen(hull):
        hull._npoints = 0
        hull.mins = np.full(hull.sample.nvar, np.inf)
        hull.maxs = np.full(hull.sample.nvar, -np.inf)
        
        hull.bp = np.full(len(hull.direct_plan), -np.inf)
        hull.bm = np.full(len(hull.direct_plan), np.inf)
        
        #hull.update()
    
    #č nejsem jist, jestli ten update vůbec dělat.
    #č lze navrhnout třidu aj tak, že sama bude hlídat změny.
    #č Jenomže, co kdybychom ten automatický update nechtěli?
    def _update(hull):
        if hull._npoints < hull.npoints:
            #hull.points = getattr(hull.sample, hull.space)
            new_points = hull.points[hull._npoints:]
            
            new_mins = np.nanmin(new_points, axis=0)
            new_maxs = np.nanmax(new_points, axis=0)
            
            hull.mins = np.nanmin((new_mins, hull.mins), axis=0)
            hull.maxs = np.nanmax((new_maxs, hull.maxs), axis=0)
            
            A = hull.direct_plan @ new_points.T
            new_bp = np.nanmax(A, axis=1)
            new_bm = np.nanmin(A, axis=1)
            
            hull.bp = np.nanmax((new_bp, hull.bp), axis=0)
            hull.bm = np.nanmin((new_bm, hull.bm), axis=0)
            
            hull._npoints = hull.npoints

    @property
    def points(hull):
        return getattr(hull.sample, hull.space)
        
    @property
    def npoints(hull):
        return len(hull.sample)
        
    @property
    def nsimplex(hull):
        return 2 * (len(hull.direct_plan) + hull.sample.nvar)
        
    @property
    def A(hull):
        hull._update()
        #č žádná optimizace, ale kdo tu funkci bude spouštět?
        diag = np.diag(np.ones(hull.sample.nvar))
        return np.vstack((diag, -diag, -hull.direct_plan, hull.direct_plan))
        
    @property
    def b(hull):
        hull._update()
        return np.concatenate((-hull.maxs, hull.mins, hull.bm, -hull.bp))
        
    @property
    def equations(hull):
        hull._update()
        #č žádná optimizace, ale kdo tu funkci bude spouštět?
        diag = np.diag(np.ones(hull.sample.nvar))
        
        A = np.vstack((diag, -diag, -hull.direct_plan, hull.direct_plan))
        b = np.concatenate((-hull.maxs, hull.mins, hull.bm, -hull.bp))
        return np.hstack((A,b[:,None]))

    def is_inside(hull, nodes):
        hull._update()
        x = getattr(nodes, hull.space)
        
        more = np.all(x < hull.maxs, axis=1)
        less = np.all(x > hull.mins, axis=1)
        
        
        #E [normal, offset] forming the hyperplane equation of the facet (see Qhull documentation for more)
        A = hull.direct_plan
        bp = np.atleast_2d(hull.bp).T
        bm = np.atleast_2d(hull.bm).T
        
        # N=ns, E - number of hyperplane equations
        ExN = A @ x.T
        higher = np.all(ExN < bp, axis=0)
        lower = np.all(ExN > bm, axis=0)
        return np.all((more, less, higher, lower), axis=0)
    
    def is_outside(hull, nodes): 
        return ~hull.is_inside(nodes)
        
    def get_design_points(hull):
        hull._update()
        sample_model = -hull.A * hull.b.reshape(-1,1)
        return hull.sample.f_model.new_sample(sample_model, space=hull.space)
        
    # shortcut for Ghull
    # valid only if space==G
    def get_r(hull):
        if hull.space=='G':
            hull._update()
            return min(-np.max(hull.mins), np.min(hull.maxs), -np.max(hull.bm), np.min(hull.bp))
        else:
            return 0
    
            

class QHull:
    # take some global function
    _orth_helper = _orth_helper
    #get_orth_outside = get_orth_outside
    get_orth_equations = get_orth_equations
    
    _2FORM_helper = _2FORM_helper
    #get_2FORM_outside = get_2FORM_outside
    get_2FORM_equations = get_2FORM_equations
    
    def __init__(self, sample, space='G', incremental=True):
        self.sample = sample
        self.incremental = incremental
        self.space = space
    
    def regen(self):
        points = getattr(self.sample, self.space)
        self.convex_hull = spatial.ConvexHull(points, incremental=self.incremental)
        
    
    def __getattr__(self, attr):
        #č branime se rekurzii
        # defend against recursion
        #оӵ рекурсилы пезьдэт!
        if attr == 'convex_hull':
            #č zkusme rychle konvexní obálky sestavit
            #č a ihned ji vrátit
            self.regen()
            return self.convex_hull
        
        elif attr == 'enough_points':
            return self.sample.nvar < self.sample.nsim
            
        elif attr == 'A':
            return self.convex_hull.equations[:,:-1]
        elif attr == 'b':
            return self.convex_hull.equations[:,-1]
            
        elif attr == 'points':
            if self.enough_points:
                return self.convex_hull.points
            else:
                return getattr(self.sample, self.space)
                
        elif attr == 'npoints':
            if self.enough_points:
                return self.convex_hull.npoints
            else:
                return len(self.sample)
        
        elif attr == 'nsimplex':
            if self.enough_points:
                return self.convex_hull.nsimplex
            else:
                return 0
            
            
        #ё По всем вопросам обращайтесь 
        #ё на нашу горячую линию    
        else:
            self._update() #č dycky čerstý chleba!
            return getattr(self.convex_hull, attr)
            
    
    #č nejsem jist, jestli ten update vůbec dělat.
    #č lze navrhnout třidu aj tak, že sama bude hlídat změny.
    #č Jenomže, co kdybychom ten automatický update nechtěli?
    def _update(self):
        if self.convex_hull.npoints < self.sample.nsim:
            if self.incremental:
                points = getattr(self.sample, self.space)  
                self.convex_hull.add_points(points[self.convex_hull.npoints:])
            else:
                self.regen()
        
        
    def is_inside(self, nodes):
        if self.enough_points:
            self._update()
            x = getattr(nodes, self.space)
            
            #E [normal, offset] forming the hyperplane equation of the facet (see Qhull documentation for more)
            A = self.convex_hull.equations[:,:-1]
            b = self.convex_hull.equations[:,-1]
            
            # N=ns, E - number of hyperplane equations
            ExN = A @ x.T + np.atleast_2d(b).T
            mask = np.all(ExN < 0, axis=0)
            return mask
        else:
            return np.full(len(nodes), False)
    
    def is_outside(hull, nodes): 
        return ~hull.is_inside(nodes)
        
    def get_design_points(hull):
        hull._update()
        sample_model = -hull.A * hull.b.reshape(-1,1)
        return hull.sample.f_model.new_sample(sample_model, space=hull.space)
        
    def get_design_point(hull):
        hull._update()
        to_fire = np.nanargmax(hull.b) #č tak, pro jistotu
        sample_model = -hull.A[to_fire] * hull.b[to_fire]
        return hull.sample.f_model.new_sample(sample_model, space=hull.space)
    
        
    # shortcut for Ghull
    # valid only if space==G
    def get_r(hull):
        if hull.space=='G':
            if hull.enough_points:
                hull._update()
                b = hull.convex_hull.equations[:,-1]
                return -np.nanmax(b)
            else:
                return 0
        else:
            return 0
    
    def shot(hull, ns, use_MC=False):
        try:
            # take global function for shot()
            return shot(hull, ns, use_MC)
        except:
            pass
            
    def fire(hull, ns, use_MC=False):
        try:
            # take global function for fire()
            return fire(hull, ns, use_MC)
        except:
            pass
    
    def get_orth_outside(hull):
        try:
            return get_orth_outside(hull)
        except:
            #č asi máme problém s triangulací 
            #č odvažně vratíme jedničku
            #č neřikejte mi, že musím pěčlivějc zpracovavat chyby!
            return 1
    
    def get_2FORM_outside(hull):
        try:
            return get_2FORM_outside(hull)
        except:
            #č asi máme problém s triangulací 
            #č odvažně vratíme jedničku
            #č neřikejte mi, že musím pěčlivějc zpracovavat chyby!
            return 1
                




Mode Type Size Ref File
100644 blob 28117 0907e38499eeca10471c7d104d4b4db30b8b7084 IS_stat.py
100644 blob 6 0916b75b752887809bac2330f3de246c42c245cd __init__.py
100644 blob 73368 3d245b8568158ac63c80fa0847631776a140db0f blackbox.py
100644 blob 11243 10c424c2ce5e8cdd0da97a5aba74c54d1ca71e0d candybox.py
100644 blob 26963 5106dbfd64b5403599d51598c03d9fad60956286 convex_hull.py
100644 blob 84851 811ef4d818e55cbe2d0305da62fcc2e1cada359a dicebox.py
100644 blob 36930 a775d1114bc205bbd1da0a10879297283cca0d4c estimation.py
100644 blob 34394 3f0ab9294a9352a071de18553aa687c2a9e6917a f_models.py
100644 blob 31540 a577087003a885ca7499d1ee9451e703fa9d2d36 g_models.py
100644 blob 20557 521db92b6961df7a46594d147bf4c4d49e828066 ghull.py
100644 blob 42820 1092b3b9f05b11d0c53b3aa63df2460ec355085d gl_plot.py
100644 blob 2718 5d721d117448dbb96c554ea8f0e4651ffe9ac457 gp_plot.py
100644 blob 29393 96162a5d181b8307507ba2f44bafe984aa939163 lukiskon.py
100644 blob 12867 5bc331e3ef3a8a6bfeb7ebe5e1ae4f4c58e8cc06 mart.py
100644 blob 7983 75455aa723db8bab291dcf941b92b9ffdba3aef1 mart3d.py
100644 blob 1364 944337088f060aabe029b611c27145bb28de4b38 maxes.py
100644 blob 5356 faac09f784e48599ff9a67e607a8e8a990b05d80 mgraph.py
100644 blob 2004 6ea8dc8f50a656c48f786d5a00bd6398276c9741 misc.py
100644 blob 3336 ce5171f10153dc0151d824369161d65fcbbda83b mplot.py
100644 blob 1450 4849f178b588e252b8c7f6a940d2d82ad35f6914 plot.py
100644 blob 2807 1feb1d43e90e027f35bbd0a6730ab18501cef63a plotly_plot.py
100644 blob 142508 2dd1eb8becf89c34b1243d424b9cae855b128e12 qt_plot.py
100644 blob 8206 5981023118262109fca8309d9b313b521a25f88f reader.py
100644 blob 4284 a0e0b4e593204ff6254f23a67652804db07800a6 samplebox.py
100644 blob 6558 df0e88ea13c95cd1463a8ba1391e27766b95c3a5 sball.py
100644 blob 5553 bac994ae58f1df80c7f8b3f33955af5402f5a4f3 sball_old.py
100644 blob 2605 0034d2e3f14c056541888235e59127e8f28b131d schemes.py
100644 blob 21623 281aef80556b8d22842b8659f6f0b7dab0ad71af shapeshare.py
100644 blob 48537 10f90c5614e9a04f0cd9f78e75f0db4a6becb3e4 simplex.py
100644 blob 13090 2b9681eed730ecfadc6c61b234d2fb19db95d87d spring.py
100644 blob 10940 6965eabdb5599bb22773e7fef1178f9b2bb51efe stm_df.py
100644 blob 3433 3063a1b6a132cbb5440ab95f1b6af1f1ff4266ac testcases_2D.py
100644 blob 2465 d829bff1dd721bdb8bbbed9a53db73efac471dac welford.py
100644 blob 22048 4a6014ca5255aa96059ff9ed5a7e29df98d26ffc whitebox.py
Hints:
Before first commit, do not forget to setup your git environment:
git config --global user.name "your_name_here"
git config --global user.email "your@email_here"

Clone this repository using HTTP(S):
git clone https://rocketgit.com/user/iam-git/WellMet

Clone this repository using ssh (do not forget to upload a key first):
git clone ssh://rocketgit@ssh.rocketgit.com/user/iam-git/WellMet

Clone this repository using git:
git clone git://git.rocketgit.com/user/iam-git/WellMet

You are allowed to anonymously push to this repository.
This means that your pushed commits will automatically be transformed into a merge request:
... clone the repository ...
... make some changes and some commits ...
git push origin main