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 (01bb9017fc4270c524819c48ee32b38d8b93c85d) (32056 bytes) (mode 100644) (type blob)

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

import numpy as np
#from . import IS_stat
from . import sball
#from . import f_models
from scipy import stats
from scipy import spatial
from .candybox import CandyBox
from .IS_stat import PushAndPull # for Shell_IS


#č 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


#ё рельса
#č nepodařílo se mi nějak rozumně zobecnit pro libovolný prostor
#č takže Gauss
#č (každá třida implementuje zvlášť)
#def fire(design_point_G, ns):
#    to_fire_G = design_point.G[0]
#    r = np.sqrt(np.sum(np.square(to_fire_G)))
#    a = to_fire_G / r
#    fire_from = stats.norm.cdf(r)
#    t = np.linspace(fire_from, 1, ns)
#    t = stats.norm.ppf(t)
#    fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
#    return design_point.

#č 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 fire(hull, ns):
        pass
        
     #č radší nic, než tohle    
#    def fire(hull, ns): # rail
#        R2_hull = np.sum(np.square(hull.sample.G), axis=1)
#        design_index = np.nanargmax(R2_hull)
#        r = np.sqrt(np.nanmax(R2_hull))
#        a = hull.sample.G[design_index] / r
#        fire_from = stats.norm.cdf(r)
#        t = np.linspace(fire_from, 1, ns)
#        t = stats.norm.ppf(t)
#        fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
#        
#        return hull.sample.f_model.new_sample(fire_G, space='G')
        


#č 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
        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 fire(hull, ns): # 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:
    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
    
    def fire(hull, ns):
        if hull.space == 'G':
            A = hull.equations[:,:-1]
            b = hull.equations[:,-1]
            
            to_fire = np.nanargmax(b)
            a = A[to_fire]
            fire_from = stats.norm.cdf(hull.get_r())
            t = np.linspace(fire_from, 1, ns, endpoint=False)
            t = stats.norm.ppf(t)
            fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
            
            return hull.sample.f_model.new_sample(fire_G, space='G')


class CompleteHull:
    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
    
    def fire(hull, ns):
        if hull.space == 'G':
            A = hull.equations[:,:-1]
            b = hull.equations[:,-1]
            
            to_fire = np.nanargmax(b)
            a = A[to_fire]
            fire_from = stats.norm.cdf(hull.get_r())
            t = np.linspace(fire_from, 1, ns, endpoint=False)
            t = stats.norm.ppf(t)
            fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
            
            return hull.sample.f_model.new_sample(fire_G, space='G')
            

class QHull:
    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 fire(hull, ns):
        if hull.space == 'G':
            try:
                A = hull.equations[:,:-1]
                b = hull.equations[:,-1]
                
                to_fire = np.nanargmax(b) #č tak, pro jistotu
                a = A[to_fire]
                fire_from = stats.norm.cdf(hull.get_r())
                t = np.linspace(fire_from, 1, ns, endpoint=False)
                t = stats.norm.ppf(t)
                fire_G = t.reshape(-1,1) @ a.reshape(1,-1)
                
                return hull.sample.f_model.new_sample(fire_G, space='G')
            except:
                pass


#č mým úkolem při návrhu této třidy je pořádně všecko zkomplikovat.
#č Dostávám za to peníze. 
#č Takže. Udělám 3 druhů estimátorů
# convex_hull_estimation  -2: inside, -1: outside
# shell_estimation  -22: inner, -3: shell, -11: outer
# ghull_estimation  -22: inner, -21: shell inside, -12: shell outside, -11: outer 

try: # try to outthink f****** OS's memory management
    import os
    mem_bytes = os.sysconf('SC_PAGE_SIZE') * os.sysconf('SC_PHYS_PAGES')
except BaseException as e:
    mem_GB = 16
    mem_bytes = mem_GB * 1024**3 # hello, Windows!
    print("convex_hull failed to get amount of RAM installed. %s GB is supposed."% mem_GB, repr(e))
    print("BTW, you are using Windows, aren't you?")
    

class Ghull:
    def __init__(self, hull, use_MC=False, non_Gaussian_reduction=0.9):
        self.hull = hull
        self.shell = sball.Shell(hull.sample.nvar)
        self.outside_dist = sball.Shell(hull.sample.nvar)
        self.sample = hull.sample
        
        if use_MC:
            self.gint = Shell_MC(hull, self.shell, non_Gaussian_reduction)
        else:
            self.gint = Shell_IS(hull, self.shell, non_Gaussian_reduction)
        
        
    def fire(self, ns):
        nodes = self.hull.fire(ns)
        if nodes is not None:
            return nodes
        else:
            return self.boom(ns)
            
    def boom(self, ns):
        self.outside_dist.set_bounds(self.get_R())
        nodes_G = self.outside_dist.rvs(ns)
        nodes = self.sample.f_model.new_sample(nodes_G, space='G')
        return nodes
        
    def get_R(self):
        sum_squared = np.sum(np.square(self.sample.G), axis=1)
        #index = np.argmax(sum_squared)
        return np.sqrt(np.nanmax(sum_squared))
        
    def get_r(self):
        "calculates minimal signed distance to planes. Can therefore be negative"
        #return -np.nanmax(self.hull.b)
        if self.hull.space == 'G':
            return self.hull.get_r()
        else:
            #č get_design_points() nemusí být.
            #č QHull může nemít dostatek teček
            try:
                hull_points = self.hull.get_design_points()
                sum_squared = np.sum(np.square(hull_points.G), axis=1)
                hull_r = np.sqrt(np.nanmin(sum_squared))
            except BaseException as e:
                msg = "cannot get hull points from the convex hull"
                print(self.__class__.__name__ +":", msg, repr(e))
                hull_r = np.inf
                
            
            # ask our experimentation team
            gint_r = self.gint.get_r()
            
            #č podle mně, nemůže se stat, že by metoda vratila:
            #č 1. nenůlový poloměr a že by dokonce i centralní bod byl venku.
            #č 2. poloměr větší jak R
            #č Odpovědnost za to kladu na gint.
            return min(hull_r, gint_r)
     
    def get_shell_estimation(self):
        shell = self.shell
        r = self.get_r()
        R = self.get_R()
        
        if r<0:
            shell.set_bounds(0, R)
        else:
            shell.set_bounds(r, R)
        
        # shell_estimation  -22: inner, -3: shell, -11: outer
        shell_estimation = {-22:shell.ps, -3: shell.p_shell, -11: shell.pf}
        global_stats = {"nsim":self.sample.nsim, "ndim":self.sample.nvar, \
                        "nfacets": self.hull.nsimplex, "r":r, "R":R, \
                    "inner":shell.ps, "shell":shell.p_shell, "outer":shell.pf}
        return shell_estimation, global_stats
        
    #č nie    
    ##č pro mně je důležité, aby bylo možné rychle přepinat mezi metodama,
    ##č aby bylo možné výsledky jednoduše porovnavat.
    ##č Proto explicitně posílám priznak, jakou metodu použit.
    def integrate(self, nis, callback_all=None, callback_outside=None):
        #č no to teda disajn třidy je doopravdy hroznej
        # .get_shell_estimation() will calculate radiuses and will update shell
        shell_estimation, global_stats = self.get_shell_estimation()
        # shell_estimation  -22: inner, -3: shell, -11: outer
        shell_estimation.pop(-3)
        # ghull_estimation  -22: inner, -21: shell inside, -12: shell outside, -11: outer 
        ghull_estimation = shell_estimation; del(shell_estimation)
        
        
        shell_ps, shell_pf, stats = self.gint.integrate(nis, callback_all, callback_outside)
        ghull_estimation[-21] = shell_ps
        ghull_estimation[-12] = shell_pf
        global_stats.update(stats)
        
        # convex_hull_estimation  -2: inside, -1: outside
        inside = shell_ps + self.shell.ps
        outside = shell_pf + self.shell.pf
        convex_hull_estimation = {-2: inside, -1: outside}
        global_stats["inside"] = inside
        global_stats["outside"] = outside
        
        return ghull_estimation, convex_hull_estimation, global_stats
        
        
        

#č pomocná třída pro integrování
#č nejsem tedy jist, jestli je to nejlepší napad - 
#č dělit Ghull na podtřídy, delegovat funkcionalitu
#č ale jinak Ghull se stavá hodně překomplikovaným.
#č nic lepšího mně nenapadá, přemyšlel jsem dlouho.
class Shell_MC:
    def __init__(self, hull, shell, non_Gaussian_reduction=0.98):
        self.shell = shell
        self.hull = hull
        self.nvar = hull.sample.nvar
        
        self.integration_cutoff = np.inf
        self.non_Gaussian_reduction = non_Gaussian_reduction
        
        # test if r>0
        self.r = 0
        central_G = np.full(hull.sample.nvar, 0, dtype=np.int8)
        self.DP = self.hull.sample.f_model.new_sample(central_G, space='G')
        
    def reset(self): # clear
        self.r_needed = (self.hull.space!='G')
        
        self.nsampled = 0
        self.nfailed = 0
        
    @property
    def DP_is_valid(self):
        #č dalo by se kontrolovat změny obálky 
        #č a provadět kontrolu DP pouze po změně obálky.
        #č Nejsem liný, jeden bodík prostě za to doopravdy nestoji.
        mask = self.hull.is_outside(self.DP)
        return np.any(mask)
        
    #č na rozdil od rodičovské třídy
    #č vrací odhad r na základě předchozích integrací
    #č metoda je navržena tak, aby Shell_IS jú mohl zdědit. 
    def get_r(self):
        #č bojim sa, rerukce bude aplikována dycky
        #č Bacha, metoda bude vracet nuly pro obálky v Gaussovskem prostoru!
        return self.r * self.non_Gaussian_reduction
        
    # stateless
    def rvs(self, nsampled, seats, ns): 
        "Generování vzorků (kandidátů a integračních bodů)"
        # rand_dir: prepare ns random directions on a unit d-sphere
        rand_dir = sball.get_random_directions(seats, self.nvar) #random directions
        
        # generate sampling probabilites
        left = (ns-nsampled) / ns
        right = (ns-nsampled-seats) / ns
        #č přidáme trochu zmatku.
        #č globálně jdeme r->R, localně ale R_i+ -> R_i-
        #č menší poloměry musejí jít dřív - na to zavazano nonGaussian _r
        #č převracení lokalně umožní vždycky mít alespoň jeden bodík outside,
        #č což je taky velmi vhodné vzhledem k tomu, že se _r bere z outside bodíků
        p = np.linspace(right, left, seats, endpoint=False) # probabilities for the radius
        
        # convert probabilitites into random radii
        # (distances from origin that are greater than r and less than R)
        r = self.shell.isf(p) # actually, it is the same as CDF inverse
        
        #finally a random sample from the optimal IS density:
        sample_G = rand_dir*r[:,None]
        
        return sample_G
    
    # bus analogy
    def _process_part(self, seats, nis, callback_all=None, callback_outside=None):
        # boarding
        nodes_G = self.rvs(self.nsampled, seats, nis)
        nodes = self.hull.sample.f_model.new_sample(nodes_G, space='G')
        # who is back?
        mask = self.hull.is_outside(nodes)
        if self.r_needed and np.any(mask):
            #č rvs má vzorkovat od měnšího poloměru k většímu.
            #č to znamená, že první outside bude mít nejměnší r vůbec
            self._r_check_out(nodes[mask])
            self.r_needed = False
        
        if callback_all is not None:
            # -2 = 'inside' -1 = 'outside'
            candy_nodes = CandyBox(nodes, event_id=mask-2, is_outside=mask)
            callback_all(candy_nodes)
            
        if (callback_outside is not None) and np.any(mask):
            callback_outside(nodes[mask])
        
        assert len(mask) == seats
        #č nevím proč, ale kdysi mě to vyšlo rychlejc jak obyčejný součet
        self.nfailed += len(mask[mask]) 
        self.nsampled += len(mask) 
        
    #č metoda je navržena tak, aby Shell_IS jú mohl zdědit. 
    def _r_check_out(self, outside_nodes):
        sum_squares = np.sum(np.square(outside_nodes.G), axis=1)
        arg = np.nanargmin(sum_squares)
        new_r = np.sqrt(sum_squares[arg])
        if (not self.DP_is_valid) or (new_r < self.r):
            self.DP = outside_nodes[arg]
            self.r = new_r
        
    
    #č metoda je navržena tak, aby Shell_IS jú mohl zdědit. 
    def integrate(self, nis, callback_all=None, callback_outside=None):
        self.reset()
        
        if self.hull.nsimplex == 0:
            bus = self.integration_cutoff
        else:
            bus = int(mem_bytes / self.hull.nsimplex / 8 / 3) + 1
        while self.nsampled < nis:
            
            seats = min(nis - self.nsampled, self.integration_cutoff, bus)
            try: 
                self._process_part(seats, nis, callback_all, callback_outside)
            except MemoryError as m:
                print(self.__class__.__name__ +":", "memory error, %s sedaček" % seats, repr(m))
                self.integration_cutoff = int(np.ceil(seats/2))
        
        assert nis == self.nsampled
        
        return self._get_result()
        
        
    def _get_result(self):
        
        nis = self.nsampled
        nf = self.nfailed
        ns = nis - nf
        p_shell = self.shell.p_shell
        shell_pf = nf/nis * p_shell
        shell_ps = ns/nis * p_shell
        
        stats = dict()
        stats["shell_budget"] = nis
        stats["shell_inside"] = shell_ps
        stats["shell_outside"] = shell_pf
        
        return shell_ps, shell_pf, stats



class Shell_IS(Shell_MC):
    def reset(self): # clear
        self.r_needed = (self.hull.space!='G')
        
        self.nsampled = 0
        self.pp = PushAndPull()
        
        
    # stateless
    def rvs(self, nsampled, seats, ns): 
        "Generování vzorků (kandidátů a integračních bodů)"
        # rand_dir: prepare ns random directions on a unit d-sphere
        rand_dir = sball.get_random_directions(seats, self.nvar) #random directions
        
        #č poloměry bereme ze skořapky
        #č za správné nastavení (stejně jako u MC)
        #č zodpovidá uživatel třídy
        #č třída vlastní odhad r nijak nevyuživá!
        r = self.shell.r
        R = self.shell.R
        delta_r = R - r
        left = nsampled / ns * delta_r + r
        right = (nsampled + seats) / ns * delta_r + r
        #č přidáme trochu zmatku.
        #č globálně jdeme r->R, localně ale R_i+ -> R_i-
        #č menší poloměry musejí jít dřív - na to zavazano nonGaussian _r
        #č převracení lokalně umožní vždycky mít alespoň jeden bodík outside,
        #č což je taky velmi vhodné vzhledem k tomu, že se _r bere z outside bodíků
        rs = np.linspace(right, left, seats, endpoint=False)
        
        #finally a random sample from the optimal IS density:
        sample_G = rand_dir*rs[:,None]
        
        #č a jsme doma. Platba za špatný design.
        #č Potřebujeme hustotu 1D rozdělení, implementovanou v Radial
        #č Shell se síce dědí od Radial, ale implementuje .pdf() jako sdruženou 
        #č hustotu v nD Gaussovskem prostoru
        #
        #č vahy jsou definovány jako podil původního rozdělení k vzorkovácímu
        #č původní - Chi rozdělení. nemá zde být žádný p_shell
        f_pdf = stats.chi.pdf(rs, self.nvar)
        #č vzorkovácí. 
        #h_pdf = 1/delta_r
        weights = f_pdf * delta_r
        
        return sample_G, weights
    
    # bus analogy
    def _process_part(self, seats, nis, callback_all=None, callback_outside=None):
        # boarding
        nodes_G, weights = self.rvs(self.nsampled, seats, nis)
        nodes = self.hull.sample.f_model.new_sample(nodes_G, space='G')
        # who is back?
        mask = self.hull.is_outside(nodes)
        if self.r_needed and np.any(mask):
            #č rvs má vzorkovat od měnšího poloměru k většímu.
            #č to znamená, že první outside bude mít nejměnší r vůbec
            self._r_check_out(nodes[mask])
            self.r_needed = False
        
        assert len(mask) == seats
        self.pp.add(weights, mask)
        self.nsampled += len(mask) 
        
        if callback_all is not None:
            # -2 = 'inside' -1 = 'outside'
            candy_nodes = CandyBox(nodes, event_id=mask-2, is_outside=mask, weights=weights)
            callback_all(candy_nodes)
            
        if (callback_outside is not None) and np.any(mask):
            callback_outside(nodes[mask])
        
        
    def _get_result(self):
        
        nis = self.nsampled
        #č nejdřív true, pak false. Posilali jsme is_outside
        mean_pf, mean_ps = self.pp.mean
        var_pf, var_ps = self.pp.var
        shell_pf, shell_ps = self.pp.correct_means(p_overall=self.shell.p_shell)
        
        stats = dict()
        stats["shell_budget"] = nis
        stats["shell_inside_measured"] = mean_ps
        stats["shell_outside_measured"] = mean_pf
        stats["shell_inside_var"] = var_ps
        stats["shell_outside_var"] = var_pf
        stats["shell_inside"] = shell_ps
        stats["shell_outside"] = shell_pf
        
        return shell_ps, shell_pf, stats


Mode Type Size Ref File
100644 blob 26998 eb786c41e6db077a9fd4999e8bba97eed20d2f50 IS_stat.py
100644 blob 6 0916b75b752887809bac2330f3de246c42c245cd __init__.py
100644 blob 73368 3d245b8568158ac63c80fa0847631776a140db0f blackbox.py
100644 blob 11243 10c424c2ce5e8cdd0da97a5aba74c54d1ca71e0d candybox.py
100644 blob 32056 01bb9017fc4270c524819c48ee32b38d8b93c85d convex_hull.py
100644 blob 84177 a7b9d900b153c4f6339281571fe4b32d7ead92df dicebox.py
100644 blob 36930 a775d1114bc205bbd1da0a10879297283cca0d4c estimation.py
100644 blob 34394 3f0ab9294a9352a071de18553aa687c2a9e6917a f_models.py
100644 blob 31540 a577087003a885ca7499d1ee9451e703fa9d2d36 g_models.py
100644 blob 42820 1092b3b9f05b11d0c53b3aa63df2460ec355085d gl_plot.py
100644 blob 2718 5d721d117448dbb96c554ea8f0e4651ffe9ac457 gp_plot.py
100644 blob 29393 96162a5d181b8307507ba2f44bafe984aa939163 lukiskon.py
100644 blob 12028 dabcfd7eb6c467ff25efa47eccebfd21697c9473 mart.py
100644 blob 7983 75455aa723db8bab291dcf941b92b9ffdba3aef1 mart3d.py
100644 blob 5356 faac09f784e48599ff9a67e607a8e8a990b05d80 mgraph.py
100644 blob 2004 6ea8dc8f50a656c48f786d5a00bd6398276c9741 misc.py
100644 blob 2439 fe482f41cb07d6d8a079553aa09b13a8a82d512d mplot.py
100644 blob 1450 4849f178b588e252b8c7f6a940d2d82ad35f6914 plot.py
100644 blob 2807 1feb1d43e90e027f35bbd0a6730ab18501cef63a plotly_plot.py
100644 blob 138773 fabb2365e733b98697b5e3ce8943db9443747c5f qt_plot.py
100644 blob 8206 5981023118262109fca8309d9b313b521a25f88f reader.py
100644 blob 4284 a0e0b4e593204ff6254f23a67652804db07800a6 samplebox.py
100644 blob 6396 8dee4de0d2afced5c7a5c3e91684a9e139b5aa00 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