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".

/gl_plot.py (cbecebc34eb21cc476b640b6c95ca9ba5ce745e1) (52315 bytes) (mode 100644) (type blob)

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

import pyqtgraph as pg
import pyqtgraph.opengl as gl
from pyqtgraph.Qt import QtGui
from pyqtgraph.Qt import QtCore


import numpy as np
import pandas as pd # required for estimation graph

from scipy import spatial

from . import estimation as stm
from . import qt_plot
         
def qt_gui_plot_3d(sample_box, space='R', *args, **kwargs):
    """
    This function will start Qt graph window in event loop
    """
    
    return QtGuiPlot3D(sample_box, space, *args, **kwargs)
    

### Define a top-level widget to hold everything
class QtGuiPlot3D(qt_plot.QtGuiPlot2D):
    #sb_runned = QtCore.pyqtSignal()
    
    def initialize_central_widget(self):
        self.central_widget = gl.GLViewWidget()
        self.redraw_called.connect(self.redraw_plot)
        
        self.base_plotting = BasePlotting(self)
        self.px_size = 5
        try:
            self.space = self.sample_box.tri_space
        except AttributeError:
            pass
            
        #self.plot_widget_2d()

    def redraw_plot(self):
        try:
            self.central_widget.clear()
            self.central_widget.reset()
        except: #č 0.10 ještě neměla ty funkce
            for item in self.central_widget.items:
                item._setView(None)
            self.central_widget.items = []
            self.central_widget.update()
        self.central_widget.opts['fov'] = 3

        

    def plot_setup(self):
        self.view_items = []
        self.view_items.append(UnitCube(self))
        self.view_items.append(Axes(self))
        self.view_items.append(Grid(self))
        self.view_items.append(LastShot(self))
        self.view_items.append(Density(self))
        self.view_items.append(ConvexHull(self))
        self.view_items.append(Facets(self))
        self.view_items.append(Wireframe(self))
        

        dock = dock_r = QtGui.QDockWidget("Simplex-based pf estimation", self)
        dock.setWidget(SimplexEstimationWidget(self, dock))
        self.view.addAction(dock.toggleViewAction())
        self.dockables.append(dock)
        self.addDockWidget(QtCore.Qt.RightDockWidgetArea, dock)
        
        #!
        dock = QtGui.QDockWidget("Tesselation-based pf estimation", self)
        dock.setWidget(VoronoiEstimationWidget(self, dock))
        self.view.addAction(dock.toggleViewAction())
        self.dockables.append(dock)
        #self.addDockWidget(QtCore.Qt.RightDockWidgetArea, dock)
        self.tabifyDockWidget(dock_r, dock)
        
        
        dock = QtGui.QDockWidget("Blackbox's candidates", self)
        dock.setWidget(CandidatesWidget(self, dock))
        self.view.addAction(dock.toggleViewAction())
        self.dockables.append(dock)
        self.tabifyDockWidget(dock_r, dock)


        
        



class BasePlotting(gl.GLScatterPlotItem):
    def __init__(self, w):
        super().__init__()
        self.w = w
        #self.w.box_runned.connect(self.redraw) #č dublikuje slice_changed
        self.w.space_changed.connect(self.plot)
        self.w.slice_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)
        
    #č důležité je, že GLScatterPlotItem nemá .redraw()
    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)
        
    def plot(self):
        nsim = self.w.slider.value()
        
        pos = getattr(self.w.sample_box, self.w.space)[:nsim]
        color = np.empty((nsim, 4))
        failsi = self.w.sample_box.failsi[:nsim]
        color[failsi] = (1, 0, 0, 1)
        color[~failsi] = (0, 1, 0, 1)

        self.setData(pos=pos, color=color, size=self.w.px_size*1.5)





class UnitCube(gl.GLBoxItem):
    def __init__(self, w):
        super().__init__()
        self.w = w
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.item = QtGui.QListWidgetItem('Unit cube')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.item)
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)

    def plot(self):
        if self.item.checkState() and (self.w.space in ('P', 'U')):
            self.setSize(1, 1, 1)
            self.show()
        elif self.item.checkState() and (self.w.space in ('aP', 'aU')):
            x, y, z, *__ = (*self.w.sample_box.alpha,)
            self.setSize(x, y, z)
            self.show()
        else:
            self.hide()


class LastShot(gl.GLScatterPlotItem):
    def __init__(self, w):
        super().__init__()
        self.w = w
        self.w.box_runned.connect(self.plot)
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.item = QtGui.QListWidgetItem('Last shot')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.item)
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)

    def plot(self):
        if self.item.checkState() and (self.w.last_shot is not None):
            pos = getattr(self.w.last_shot, self.w.space)
            self.setData(pos=pos, color=(0, 0, 1, .7), size=self.w.px_size*2)
            self.show()
        else:
            self.hide()



#č Nekreslí to co má (vidím pouze jednu čáru)

##class Axes(gl.GLAxisItem):
##    def __init__(self, w):
##        super().__init__()
##        self.w = w
##        self.w.redraw_called.connect(self.redraw)
##
##        self.item = QtGui.QListWidgetItem('Axes')
##        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
##        self.item.setCheckState(QtCore.Qt.Checked)
##        self.w.list_view.addItem(self.item)
##        self.w.list_view.itemChanged.connect(self.plot)
##        
##
##    def redraw(self):
##        self.plot()
##        self.w.central_widget.addItem(self)
##
##    def plot(self):
##        if self.item.checkState():
##            self.show()
##        else:
##            self.hide()


class Axes:
    def __init__(self, w):
        pos = np.array([[0,0,0],[1,0,0], [0,0,0],[0,1,0], [0,0,0],[0,0,1]])
        color = np.array([[255, 0, 0, 255],[0, 0, 0, 0], [0, 255, 0, 255],[0, 0, 0, 0], [0, 0, 255, 255],[0, 0, 0, 0]])
        self.axes = gl.GLLinePlotItem(pos=pos, color=color, width=2, mode='lines')
        
        self.w = w
        self.w.redraw_called.connect(self.redraw)

        self.item = QtGui.QListWidgetItem('Axes')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.item)
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self.axes)

    def plot(self):
        if self.item.checkState():
            self.axes.show()
        else:
            self.axes.hide()



class Grid:
    def __init__(self, w):
        self.gx = gl.GLGridItem()
        self.gx.setSize(10, 10, 1)
        self.gx.rotate(90, 0, 1, 0)
        #gx.translate(-5, 0, 5)
        
        self.gy = gl.GLGridItem()
        self.gy.setSize(10, 10, 1)
        self.gy.rotate(90, 1, 0, 0)
        #gy.translate(0, -10, 10)
        
        self.gz = gl.GLGridItem()
        self.gz.setSize(10, 10, 1)
        #self.gz.translate(0, 0, 0)
        
        self.w = w
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.item = QtGui.QListWidgetItem('Grid')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Unchecked)
        self.w.list_view.addItem(self.item)
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self.gx)
        self.w.central_widget.addItem(self.gy)
        self.w.central_widget.addItem(self.gz)

    def plot(self):
        if self.item.checkState() and (self.w.space in ( 'Rn', 'GK', 'G')):
            self.gx.show()
            self.gy.show()
            self.gz.show()
        else:
            self.gx.hide()
            self.gy.hide()
            self.gz.hide()




class ConvexHull(gl.GLMeshItem):
    def __init__(self, w):
        ## Mesh item will automatically compute face normals.
        super().__init__(smooth=True, shader='edgeHilight') # edgeHilight shaded
        #self.setGLOptions('opaque')
        self.nsim = -1
        
        self.w = w
        self.w.box_runned.connect(self.plot)
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.failure_item = QtGui.QListWidgetItem('ConvexHull Failure')
        self.failure_item.setFlags(self.failure_item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.failure_item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.failure_item)

        self.success_item = QtGui.QListWidgetItem('ConvexHull Success')
        self.success_item.setFlags(self.success_item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.success_item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.success_item)

        
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)

    def recalculate(self):
        #č velmi krátko - ConvexHull rozmichá vzorky a nadělá ďupy
        sampled_plan_tri = getattr(self.w.sample_box, self.w.sample_box.tri_space)
        tree = spatial.cKDTree(sampled_plan_tri)
        self.points = self.w.sample_box.convex_hull.points
        dd, ii = tree.query(self.points, k=1, p=2)

        self.simplices = self.w.sample_box.convex_hull.simplices
        box_facets = ii[self.simplices.flatten()].reshape(-1, 3)
        self.events = self.w.sample_box.get_events(box_facets)


        nodes_colors = np.empty((self.w.sample_box.nsim, 4))
        nodes_colors[self.w.sample_box.failsi] = np.array([253, 93, 97, 0])/255 
        nodes_colors[~self.w.sample_box.failsi] = np.array([67, 255, 81, 0])/255
        

        self.vertex_colors = nodes_colors[ii]
        
        # marker
        self.nsim = self.w.sample_box.nsim


    def plot(self):
        if (self.failure_item.checkState() + self.success_item.checkState())\
                               and (self.w.space == self.w.sample_box.tri_space):
            try:
                if self.nsim != self.w.sample_box.nsim:
                    self.recalculate()
                    
                #face_colors = np.empty((len(self.simplices), 4), dtype=np.int16)
                # sorry for that
                if self.failure_item.checkState():
                    mask_1 = self.events==1
                else:
                    mask_1 = self.events==100500

                if self.success_item.checkState():
                    mask_0 = self.events==0
                else:
                    mask_0 = self.events==100500

                #mask_2 = self.events!=2

                mask = np.any([mask_0, mask_1], axis=0)

                
                self.setMeshData(vertexes=self.points, faces=self.simplices[mask], vertexColors=self.vertex_colors, \
                                 drawEdges=True, edgeColor=(1, 1, 0, 1))
            except BaseException as e:
                msg = "nepovedlo se nám spočítat konvexní obálku "
                error_msg = self.__class__.__name__ + ": " + msg + repr(e)
                print(error_msg)
                
            self.show()
        else:
            self.hide()

        

class Facets(gl.GLMeshItem):
    def __init__(self, w):
        ## Mesh item will automatically compute face normals.
        super().__init__(smooth=True, shader='shaded') # edgeHilight shaded
        #self.setGLOptions('opaque')
        self.nsim = -1
        
        self.w = w
        self.w.box_runned.connect(self.plot)
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.failure_item = QtGui.QListWidgetItem('Failure facets')
        self.failure_item.setFlags(self.failure_item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.failure_item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.failure_item)

        self.success_item = QtGui.QListWidgetItem('Success facets')
        self.success_item.setFlags(self.success_item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.success_item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.success_item)

        
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)

    def recalculate(self):
        sampled_plan_tri = getattr(self.w.sample_box, self.w.sample_box.tri_space)
        self.points = self.w.sample_box.tri.points

        self.simplices = self.w.sample_box.tri.simplices
        self.s_events = self.w.sample_box.get_events(self.simplices)

        self.mixed_simplices = self.simplices[self.s_events==2]

        self.facets = np.vstack((self.mixed_simplices[:,:3], self.mixed_simplices[:,1:],\
                                 self.mixed_simplices[:,[0,2,3]], self.mixed_simplices[:,[0,1,3]]))

        # f_events
        self.events = self.w.sample_box.get_events(self.facets)
                    

        self.nodes_colors = np.empty((self.w.sample_box.nsim, 4))
        self.nodes_colors[self.w.sample_box.failsi] = np.array([203, 83, 87, 0])/255 
        self.nodes_colors[~self.w.sample_box.failsi] = np.array([57, 205, 71, 0])/255
        
        
        # marker
        self.nsim = self.w.sample_box.nsim


    def plot(self):
        if (self.failure_item.checkState() + self.success_item.checkState())\
                               and (self.w.space == self.w.sample_box.tri_space):
            try:
                if self.nsim != self.w.sample_box.nsim:
                    self.recalculate()
                    
                #face_colors = np.empty((len(self.simplices), 4), dtype=np.int16)
                # sorry for that
                if self.failure_item.checkState():
                    mask_1 = self.events==1
                else:
                    mask_1 = self.events==100500

                if self.success_item.checkState():
                    mask_0 = self.events==0
                else:
                    mask_0 = self.events==100500

                #mask_2 = self.events!=2

                mask = np.any([mask_0, mask_1], axis=0)

                self.setMeshData(vertexes=self.points, faces=self.facets[mask], vertexColors=self.nodes_colors, \
                                 drawEdges=True, edgeColor=(1, 1, 0, 1))
            except BaseException as e:
                msg = "nepovedlo se nám spočítat stěny "
                error_msg = self.__class__.__name__ + ": " + msg + repr(e)
                print(error_msg)
                
            self.show()
        else:
            self.hide()   
        
        


class Wireframe(gl.GLMeshItem):
    def __init__(self, w):
        super().__init__(smooth=False, drawFaces=False, drawEdges=True) 
        self.nsim = -1
        
        self.w = w
        self.w.box_runned.connect(self.plot)
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        self.item = QtGui.QListWidgetItem('Wireframe')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Checked)
        self.w.list_view.addItem(self.item)

        
        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self)

    def recalculate(self):
        sampled_plan_tri = getattr(self.w.sample_box, self.w.sample_box.tri_space)
        self.points = self.w.sample_box.tri.points

        self.simplices = self.w.sample_box.tri.simplices


        self.facets = np.vstack((self.simplices[:,:3], self.simplices[:,1:],\
                                 self.simplices[:,[0,2,3]], self.simplices[:,[0,1,3]]))



        self.setMeshData(vertexes=self.points, faces=self.facets, drawEdges=True, edgeColor=(0.5, 0.5, 0.5, 1))
        
        # marker
        self.nsim = self.w.sample_box.nsim


    def plot(self):
        if self.item.checkState() and (self.w.space == self.w.sample_box.tri_space):
            try:
                if self.nsim != self.w.sample_box.nsim:
                    self.recalculate()

                
            except BaseException as e:
                msg = "nepovedlo se nám spočítat stěny "
                error_msg = self.__class__.__name__ + ": " + msg + repr(e)
                print(error_msg)
                
            self.show()
        else:
            self.hide()   
        


#☺ čota mne ne chočetsa eto dělať
            
##class FlexFrame(GLLinePlotItem):
##    def __init__(self, w):
##        super().__init__(smooth=False, drawFaces=False, drawEdges=True) 
##        self.nsim = -1
##        
##        pos = np.array([[0,0,0],[1,0,0], [0,0,0],[0,1,0], [0,0,0],[0,0,1]])
##        color = np.array([[255, 0, 0, 255],[0, 0, 0, 0], [0, 255, 0, 255],[0, 0, 0, 0], [0, 0, 255, 255],[0, 0, 0, 0]])
##        self.axes = gl.GLLinePlotItem(pos=pos, color=color, width=2, mode='lines')
##        
##        self.w = w
##        self.w.slice_changed.connect(self.plot)
##        self.w.space_changed.connect(self.plot)
##        self.w.redraw_called.connect(self.redraw)
##
##        self.item = QtGui.QListWidgetItem('FlexFrame')
##        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
##        self.item.setCheckState(QtCore.Qt.Checked)
##        self.w.list_view.addItem(self.item)
##        self.w.list_view.itemChanged.connect(self.plot)
##
##
##    def recalculate(self):
##        self.points = self.w.sample_box.tri.points
##
##        self.simplices = self.w.sample_box.tri.simplices
##
##
##        self.facets = np.vstack((self.simplices[:,:3], self.simplices[:,1:],\
##                                 self.simplices[:,[0,2,3]], self.simplices[:,[0,1,3]]))
##
##
##
##        self.setMeshData(vertexes=self.points, faces=self.facets, drawEdges=True, edgeColor=(0.5, 0.5, 0.5, 1))
##
##        
##        
##        # take coordinates in the space, where triangulation has been performed
##        sampled_plan_tri = getattr(self.sample_box, self.sample_box.tri_space)
##        ns = 100
##        
##        if len(self.triangulation) < len(self.sample_box.tri.simplices):
##            x = y = ()
##            gap = len(self.sample_box.tri.simplices) - len(self.triangulation)
##            for __ in range(gap):
##                self.triangulation.append(self.plotWidget.plot(x, y, pen=0.7))
##        
##        for triangle, plot_item in zip(self.sample_box.tri.simplices, self.triangulation):
##            x_tri = np.linspace(sampled_plan_tri[triangle[0],0], sampled_plan_tri[triangle[1],0], ns, endpoint=False)
##            y_tri = np.linspace(sampled_plan_tri[triangle[0],1], sampled_plan_tri[triangle[1],1], ns, endpoint=False)
##            x_tri = np.append(x_tri, np.linspace(sampled_plan_tri[triangle[1],0], sampled_plan_tri[triangle[2],0], ns, endpoint=False))
##            y_tri = np.append(y_tri, np.linspace(sampled_plan_tri[triangle[1],1], sampled_plan_tri[triangle[2],1], ns, endpoint=False))
##            x_tri = np.append(x_tri, np.linspace(sampled_plan_tri[triangle[2],0], sampled_plan_tri[triangle[0],0], ns, endpoint=True))
##            y_tri = np.append(y_tri, np.linspace(sampled_plan_tri[triangle[2],1], sampled_plan_tri[triangle[0],1], ns, endpoint=True))
##            
##            tri_bound_tri = np.array((x_tri, y_tri)).T
##            # vytvořme sample
##            tri_bound = self.sample_box.sampled_plan.new_sample(tri_bound_tri, space=self.sample_box.tri_space)
##            
##            # Update the data
##            #print("Суредасько", tri_bound)
##            plot_item.setData(*getattr(tri_bound, self.space).T)
##                    
##
##
##
##        
##        # marker
##        self.nsim = self.w.sample_box.nsim
##
##
##    def redraw(self):
##        self.plot()
##        self.w.central_widget.addItem(self.axes)
##
##    def plot(self):
##        if self.item.checkState():
##            try:
##                self.recalculate()
##            except BaseException as e:
##                msg = "error during triangulation drawing"
##                error_msg = self.__class__.__name__ + ": " + msg + repr(e)
##                print(error_msg)
##                
##            self.axes.show()
##        else:
##            self.axes.hide()




class Density:
    def __init__(self, w):
        self.w = w
        self.w.space_changed.connect(self.plot)
        self.w.redraw_called.connect(self.redraw)

        
        self.v = gl.GLVolumeItem(self.recalculate())
        self.space = self.w.space
        self.v.scale(0.1, .10, .10)
        self.v.translate(-5,-5,-5)

        self.item = QtGui.QListWidgetItem('Density')
        self.item.setFlags(self.item.flags() | QtCore.Qt.ItemIsUserCheckable)
        self.item.setCheckState(QtCore.Qt.Unchecked)
        self.w.list_view.addItem(self.item)

        self.w.list_view.itemChanged.connect(self.plot)
        

    def redraw(self):
        self.plot()
        self.w.central_widget.addItem(self.v)

    def recalculate(self):
        #č strašně drahý, nechce se mi to pořádně udělat
        #data = np.fromfunction(lambda i,j,k:self.w.sample_box.pdf((np.array([[i,j,k]])-[50, 50, 50])/10, self.w.space), (100,100,100))

        # ničemu nerozumím
        grid = (np.mgrid[0:100,0:100, 0:100].T.reshape(-1, 3) - [50, 50, 50])/10
        shape = (100, 100, 100)
        sample = np.array((grid[:, 2], grid[:, 1], grid[:, 0])).T
        data = self.w.sample_box.sample_pdf(sample, self.w.space)
        data = data.reshape((100, 100, 100))
        
        d2 = np.empty(data.shape + (4,), dtype=np.ubyte)
        d2[..., 0] = np.full(shape, 0)
        d2[..., 1] = np.full(shape, 200)
        d2[..., 2] = np.full(shape, 200)
        d2[..., 3] = data * (255./data.max())

        return d2
        
        


    def plot(self):
        if self.item.checkState():
            if self.w.space != self.space:
                self.v.setData(self.recalculate())
                self.space = self.w.space
            self.v.show()
        else:
            self.v.hide()   






"""
=============
График виӝет 
Grafy
Estimation graph widgets
========================
+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
"""
        



class TriEstimationGraph(pg.PlotWidget):
    def __init__(self, black_box, tri_estimation_name='TRI_current_estimations',   parent=None, *args, **kwargs):
        super().__init__(parent)
        
        self.black_box = black_box
        self.tri_estimation_name = tri_estimation_name
        self.setBackground('w')
        x = y = () # zde jen vytvoříme kostru, nakrmime daty v .redraw()
        
        #self.pen_f = self.plot(x, y, fillLevel=0, brush='r') # red
        pen = pg.mkPen(color=(255, 0, 0), width=3)#, style=QtCore.Qt.DashLine)
        self.pen_f = self.plot(x, y, pen=pen) # red
        self.pen_mix = self.plot(x, y, fillLevel=0, brush=(255, 165, 0)) # orange
        self.pen_outside = self.plot(x, y, fillLevel=0, brush=0.8) # grey
        self.pen_success = self.plot(x, y, fillLevel=0, brush='g') # green
        self.pen_exact = self.plot(x, y, pen='b') # blue
        
        self.setLogMode(False, True)
        
        self.redraw()
        
   
    def redraw(self):
        try: # тут всё что угодно может пойти не так
            data = self.black_box.guessbox.estimations[self.tri_estimation_name]
            x = data[0]
            # it can be effectively done with pandas
            df = pd.DataFrame(data[1])
            # -1 = 'out', 0=success, 1=failure, 2=mix
            df.rename(columns={-1:'out', 0:'success', 1:'failure', 2:'mix'}, inplace=True)
           
            # Update the data
            # když se někde objeví nula se zapnutým LogModem - qtpygraph hned spadne a není možne ten pad zachytit
            def nonzero(data): return np.where(data > 0, data, 1) 
            self.pen_f.setData(np.array(x)[df.failure.to_numpy() > 0], df.failure.to_numpy()[df.failure.to_numpy() > 0])
            self.pen_mix.setData(x, nonzero(df.failure))
            self.pen_outside.setData(x, nonzero(df.failure+df.mix))
            self.pen_success.setData(x, nonzero(df.failure+df.mix+df.out)) # kontrolu, zda je to 1, nechame uživateli
            
            self.pen_exact.setData((min(x),max(x)), (self.black_box.pf_exact, self.black_box.pf_exact))
        
        except BaseException as e:
            print(self.__class__.__name__ + ":", repr(e))
        
        
        # pen_f.opts['logMode']
        # pen_outside.setLogMode(False, False)
        #setLogMode(False, False)
        




class SimpleSimplexEstimationGraph(pg.PlotWidget):
    def __init__(self, black_box,   parent=None, *args, **kwargs):
        super().__init__(parent)
        
        self.black_box = black_box
        
    
        self.log_x_chk = QtGui.QAction("Log X", self.plotItem.ctrlMenu)
        self.log_x_chk.setCheckable(True)
        self.log_y_chk = QtGui.QAction("Log Y", self.plotItem.ctrlMenu)
        self.log_y_chk.setCheckable(True)
        # setLogMode(self, x=None, y=None)
        self.log_x_chk.triggered.connect(lambda: self.setLogMode(x=self.log_x_chk.isChecked()))
        self.log_y_chk.triggered.connect(self.setup)
        
        self.reaction = QtGui.QAction("Redraw", self.plotItem.ctrlMenu)
        self.reaction.triggered.connect(self.redraw)
        
        #pw.plotItem.ctrlMenu.actions()
        self.plotItem.ctrlMenu.removeAction(self.plotItem.ctrlMenu.actions()[0])
        
        self.plotItem.ctrlMenu.insertAction(self.plotItem.ctrlMenu.actions()[0], self.log_y_chk)
        self.plotItem.ctrlMenu.insertAction(self.plotItem.ctrlMenu.actions()[0], self.log_x_chk)
        self.plotItem.ctrlMenu.insertAction(self.plotItem.ctrlMenu.actions()[0], self.reaction)
        
        
        self.setup()
        
        
    def setup(self, *args, **kwargs):
        self.clear()
        self.setBackground('w')
        x = y = () # zde jen vytvoříme kostru, nakrmime daty v .redraw()
        
        if self.log_y_chk.isChecked():
            self.setLogMode(y=True)
            #self.pen_f = self.plot(x, y, fillLevel=0, brush='r') # red
            pen = pg.mkPen(color=(255, 0, 0), width=3)#, style=QtCore.Qt.DashLine)
            self.pen_f = self.plot(x, y, pen=pen) # red
            self.pen_mix = self.plot(x, y, fillLevel=0, brush=(255, 165, 0)) # orange
            self.pen_outside = self.plot(x, y, fillLevel=0, brush=0.8) # grey
            self.pen_success = self.plot(x, y, fillLevel=0, brush='g') # green
            self.pen_exact = self.plot(x, y, pen='b') # blue
            
            
        else:
            self.setLogMode(y=False)
            self.pen_success = self.plot(x, y, fillLevel=0, brush='g') # green
            self.pen_outside = self.plot(x, y, fillLevel=0, brush=0.8) # grey
            self.pen_mix = self.plot(x, y, fillLevel=0, brush=(255, 165, 0)) # orange
            self.pen_f = self.plot(x, y, fillLevel=0, brush='r') # red
            self.pen_exact = self.plot(x, y, pen='b') # blue
        
        
        self.redraw()
        
   
    def redraw(self):
        xmin = np.inf
        xmax = -np.inf
        tri_estimation = dict()
        try: # тут всё что угодно может пойти не так
            # kruci, ještě navic i generovať pokažde znovu...
        
            # new-style: šecko leží dohromady a každý si z toho
            # bere co chce a jak chce
            # ne že by to bylo nějak šetrný
            # estimation je slovníkem
            for estimation in self.black_box.estimations:
                # nsim musí mäť každej odhad
                # pokud nemá - je třeba jej prostě opravit
                nsim = estimation['nsim']
                try: 
                    tri_estimation[nsim] = estimation['TRI_estimation']
                    if nsim > xmax:
                        xmax = nsim
                    if nsim < xmin:
                        xmin = nsim
                except KeyError as e:
                    pass #print(self.__class__.__name__ + ":", repr(e))
            
                
            # it can be effectively done with pandas
            df = pd.DataFrame(tuple(tri_estimation.values()))
            # -1 = 'out', 0=success, 1=failure, 2=mix
            df.rename(columns={-1:'out', 0:'success', 1:'failure', 2:'mix'}, inplace=True)
            df.insert(loc=0, column='x', value=tuple(tri_estimation.keys()), allow_duplicates=False)
            df.sort_values('x', inplace=True)
            x = df.x.to_numpy()
            
            if self.log_y_chk.isChecked():
                # Update the data
                # když se někde objeví nula se zapnutým LogModem - qtpygraph hned spadne a není možne ten pad zachytit
                def nonzero(data): return np.where(data > 0, data, 1) 
                self.pen_f.setData(x[df.failure.to_numpy() > 0], df.failure.to_numpy()[df.failure.to_numpy() > 0])
                self.pen_mix.setData(x, nonzero(df.failure))
                self.pen_outside.setData(x, nonzero(df.failure+df.mix))
                self.pen_success.setData(x, nonzero(df.failure+df.mix+df.out)) # kontrolu, zda je to 1, nechame uživateli
                
                
                
            else:
                # Update the data
                self.pen_f.setData(x, df.failure)
                self.pen_mix.setData(x, df.failure+df.mix)
                self.pen_outside.setData(x, df.failure+df.mix+df.out)
                self.pen_success.setData(x, df.failure+df.mix+df.out+df.success) # kontrolu, zda je to 1, nechame uivateli
                
           
           
            if (xmax - xmin) > 0:
                if hasattr(self.black_box, 'pf_exact'):
                    # poslední. I když spadne, tak už nikomu moc nevadí
                    self.pen_exact.setData((xmin,xmax), (self.black_box.pf_exact, self.black_box.pf_exact))
        
        except BaseException as e:
            print(self.__class__.__name__ + ":", repr(e))
        
        
        # pen_f.opts['logMode']
        # pen_outside.setLogMode(False, False)
        #setLogMode(False, False)
        #f = pg.FillBetweenItem(curves[i], curves[i+1], brushes[i])
        #win.addItem(f)
        
        







class VoronoiEstimationGraph(pg.PlotWidget):
    def __init__(self, black_box, parent=None, *args, **kwargs):
        super().__init__(parent)
        
        self.black_box = black_box
        self.setBackground('w')
        
        
        self.reaction = QtGui.QAction("Redraw", self.plotItem.ctrlMenu)
        self.reaction.triggered.connect(self.redraw)
        self.plotItem.ctrlMenu.insertAction(self.plotItem.ctrlMenu.actions()[0], self.reaction)
        
        
        
        # implicitně Y je v logaritmickem měřítku
        self.setLogMode(False, True)

        x = y = () # zde jen vytvoříme kostru, nakrmíme daty v .redraw()
        
        
        # nechapu, proč těm Itemům ríkám "propíska" 
        # propíska? Их есть у нас!
        
        self.Voronoi_2_point_upper_bound = self.plot(x, y, pen='y')
        self.Voronoi_2_point_lower_bound = self.plot(x, y, pen='y')
        
        fill_color = (255, 243, 154) # let's try xkcd: dark cream (#fff39a)
        self.fill = pg.FillBetweenItem(self.Voronoi_2_point_upper_bound, self.Voronoi_2_point_lower_bound, fill_color)
        self.addItem(self.fill)
        
        self.Voronoi_2_point_failure_rate = self.plot(x, y, pen=(195,46,212))
        self.Voronoi_2_point_pure_failure_rate = self.plot(x, y, pen='m')
        self.Voronoi_failure_rate = self.plot(x, y, pen='r')
        
        self.pen_exact = self.plot(x, y, pen='b') # blue
        self.pen_one = self.plot(x, y, pen='k') # black
        
        self.redraw()
    
   
    def redraw(self):
        # kruci, ještě navic i generovať pokažde znovu...
        metrics = {'Voronoi_2_point_upper_bound':{},\
                    'Voronoi_2_point_lower_bound':{},\
                    'Voronoi_2_point_failure_rate':{},\
                    'Voronoi_2_point_pure_failure_rate':{},\
                    'Voronoi_failure_rate':{},}
        xmin = np.inf
        xmax = -np.inf
        try: # тут всё что угодно может пойти не так
            # new-style: šecko leží dohromady a každý z toho
            # bere co chce a jak chce
            # ne že by to bylo nějak šetrný
            # estimation je slovníkem
            for estimation in self.black_box.estimations:
                # nsim musí mäť každej odhad
                # pokud nemá - je třeba jej prostě opravit
                nsim = estimation['nsim']
                
                
                for metric, metric_dict in metrics.items():
                    try: 
                        if estimation[metric] > 0:
                            metric_dict[nsim] = estimation[metric]
                            if nsim > xmax:
                                xmax = nsim
                            if nsim < xmin:
                                xmin = nsim
                    except KeyError as e:
                        pass #print(self.__class__.__name__ + ":", repr(e))
            
            for metric, metric_dict in metrics.items():
                pen = getattr(self, metric)
                # nikdo neslibil, že budou v pořadí
                x = np.sort(tuple(metric_dict.keys()))
                y = np.array(tuple(metric_dict.values()))[np.argsort(tuple(metric_dict.keys()))]
                pen.setData(x, y)
                
            if (xmax - xmin) > 0:
                self.pen_one.setData((xmin,xmax), (1, 1))
                if hasattr(self.black_box, 'pf_exact'):
                    # poslední. I když spadne, tak už nikomu moc nevadí
                    self.pen_exact.setData((xmin,xmax), (self.black_box.pf_exact, self.black_box.pf_exact))
        
        except BaseException as e:
            print(self.__class__.__name__ + ":", repr(e))
        
        
        # pen_f.opts['logMode']
        # pen_outside.setLogMode(False, False)
        #setLogMode(False, False)
        #f = pg.FillBetweenItem(curves[i], curves[i+1], brushes[i])
        #win.addItem(f)







"""
=============
Эскерон виӝет 
Widgety odhadů
Estimation widgets
===================
+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
"""

        
        
class SimplexEstimationWidget(qt_plot.SimplexEstimationWidget):
        
        
    def callback(self, estimation, simplex, nodes, cell_stats, out_nodes, *args, **kwargs):
        # stm trianguľaciju pokažde provadí znovu, proto skoro nemá cenu drbat se s její znovupoužitím
        if (simplex.nvar==3):
            ns = 30
            
            
            simplex_tri = getattr(simplex, estimation['model_space'])
            x_tri = np.linspace(simplex_tri[0,0], simplex_tri[1,0], ns, endpoint=False)
            y_tri = np.linspace(simplex_tri[0,1], simplex_tri[1,1], ns, endpoint=False)
            z_tri = np.linspace(simplex_tri[0,2], simplex_tri[1,2], ns, endpoint=False)
            
            x_tri = np.append(x_tri, np.linspace(simplex_tri[1,0], simplex_tri[2,0], ns, endpoint=False))
            y_tri = np.append(y_tri, np.linspace(simplex_tri[1,1], simplex_tri[2,1], ns, endpoint=False))
            z_tri = np.append(z_tri, np.linspace(simplex_tri[1,2], simplex_tri[2,2], ns, endpoint=False))
            
            x_tri = np.append(x_tri, np.linspace(simplex_tri[2,0], simplex_tri[3,0], ns, endpoint=True))
            y_tri = np.append(y_tri, np.linspace(simplex_tri[2,1], simplex_tri[3,1], ns, endpoint=True))
            z_tri = np.append(z_tri, np.linspace(simplex_tri[2,2], simplex_tri[3,2], ns, endpoint=True))
                
            tri_bound_tri = np.array((x_tri, y_tri, z_tri)).T
            # vytvořme sample
            tri_bound = self.sb_item.sample_box.sampled_plan.new_sample(tri_bound_tri, space=estimation['model_space'])
            
            # draw 
            pos = getattr(tri_bound, self.sb_item.space)
            plot_item = gl.GLLinePlotItem(pos=pos)
            self.sb_item.central_widget.addItem(plot_item)
            
            # uložíme data
            self.triangulation.append((tri_bound, plot_item))

            #
            #čs a eště raz
            #
            simplex_tri = getattr(simplex, estimation['model_space'])
            x_tri = np.linspace(simplex_tri[1,0], simplex_tri[3,0], ns, endpoint=False)
            y_tri = np.linspace(simplex_tri[1,1], simplex_tri[3,1], ns, endpoint=False)
            z_tri = np.linspace(simplex_tri[1,2], simplex_tri[3,2], ns, endpoint=False)
            
            x_tri = np.append(x_tri, np.linspace(simplex_tri[3,0], simplex_tri[0,0], ns, endpoint=False))
            y_tri = np.append(y_tri, np.linspace(simplex_tri[3,1], simplex_tri[0,1], ns, endpoint=False))
            z_tri = np.append(z_tri, np.linspace(simplex_tri[3,2], simplex_tri[0,2], ns, endpoint=False))
            
            x_tri = np.append(x_tri, np.linspace(simplex_tri[0,0], simplex_tri[2,0], ns, endpoint=True))
            y_tri = np.append(y_tri, np.linspace(simplex_tri[0,1], simplex_tri[2,1], ns, endpoint=True))
            z_tri = np.append(z_tri, np.linspace(simplex_tri[0,2], simplex_tri[2,2], ns, endpoint=True))
                
            tri_bound_tri = np.array((x_tri, y_tri, z_tri)).T
            # vytvořme sample
            tri_bound = self.sb_item.sample_box.sampled_plan.new_sample(tri_bound_tri, space=estimation['model_space'])
            
            # draw 
            pos = getattr(tri_bound, self.sb_item.space)
            plot_item = gl.GLLinePlotItem(pos=pos)
            self.sb_item.central_widget.addItem(plot_item)
            
            # uložíme data
            self.triangulation.append((tri_bound, plot_item))



            
        
        
        #    
        # tečičky
        #

        event = cell_stats['event']
        cell_probability = cell_stats['cell_probability']
        cm = self.param.getValues()[event][0] #č očekávám tam kolor mapu
        
        #č chcu ještě na konci prekreslit s různejma barvičkama, podle obsahu pravděpodobnosti
        # zkontrolujeme probability
        if self.max_simplices[event] < cell_probability:
            self.max_simplices[event] = cell_probability
            # a hned všecko dotyčné přebarvíme podle obsahu pravděpodobnosti
            for self_simplex in self.simplices:
                if event == self_simplex[2]['event']:
                    # zde cell_probability se rovná self.p_cell_max[event]
                    # ale bacha! Nesplet se jinde!
                    #č nechť zůstane starý nazev
                    blue_intensity = self_simplex[2]['cell_probability'] / cell_probability
                    color = cm.mapToFloat(blue_intensity)
                    #č tam prostě MUSÍ být tuple
                    self_simplex[1].setData(color=tuple(color))

        
        
        
        
        # draw tečičky
        #
        pos = getattr(nodes, self.sb_item.space)
        blue_intensity = cell_probability / self.max_simplices[event]
        color = cm.mapToFloat(blue_intensity)
        #symbolSize = np.sqrt(nodes.w / min(nodes.w)) # not bad
        size = self.param.getValues()['node (pixel) size'][0]
        #č tam prostě MUSÍ být tuple
        plot_item = gl.GLScatterPlotItem(pos=pos, size=size, color=tuple(color))
        self.sb_item.central_widget.addItem(plot_item)
        
        # uložíme data
        self.simplices.append((nodes, plot_item, cell_stats))
        
        
        
        # keep the GUI responsive :)
        self.sb_item.app.processEvents()



        
            
    def recolor(self):
        for nodes, plot_item, cell_stats in self.simplices:
            event = cell_stats['event']
            cell_probability = cell_stats['cell_probability']
            cm = self.param.getValues()[event][0] #č očekávám tam kolor mapu
        
            
            blue_intensity = cell_probability / self.max_simplices[event]
            color = cm.mapToFloat(blue_intensity)
            #symbolSize = np.sqrt(nodes.w / min(nodes.w)) # not bad
            size = self.param.getValues()['node (pixel) size'][0]
            #č tam prostě MUSÍ být tuple
            plot_item.setData(size=size, color=tuple(color))

        

    def on_space_changed(self, *args, **kwargs):
        # nejdřív triangulace
        for tri_bound, plot_item in self.triangulation:
            pos = getattr(tri_bound, self.sb_item.space)
            plot_item.setData(pos=pos)

        # teď tečičky
        for nodes, plot_item, cell_stats in self.simplices:
            pos = getattr(nodes, self.sb_item.space)
            plot_item.setData(pos=pos)
            





class VoronoiEstimationWidget(qt_plot.VoronoiEstimationWidget):

        

    def on_space_changed(self, *args, **kwargs):
        # teď tečičky
        for nodes, plot_item, cell_stats in self.cells:
            pos = getattr(nodes, self.sb_item.space)
            plot_item.setData(pos=pos)
        
    
        
        
            
        
    def node_pf_coloring(self, estimation=None, nodes=None, cell_stats=None, out_nodes=None, *args, **kwargs):
        """
        if nodes and cell_stats  provided we will add them to self.cells
        otherwise function redraw items in self.cells
        """
        
        symbol_size = self.param.getValues()['node (pixel) size'][0]
        if nodes is None:
            for nodes, plot_item, cell_stats in self.cells:
                colors = self.node_pf_colors(nodes, cell_stats)
                plot_item.setData(color=colors, size=symbol_size)
                
        
        # máme nodes, tj. jedeme poprvé        
        else:
            colors = self.node_pf_colors(nodes, cell_stats)
            pos = getattr(nodes, self.sb_item.space)
            plot_item = gl.GLScatterPlotItem(pos=pos, size=symbol_size, color=colors)
            self.sb_item.central_widget.addItem(plot_item)
            
            # uložíme data
            self.cells.append((nodes, plot_item, cell_stats))
            
            # keep the GUI responsive :)
            self.sb_item.app.processEvents()
            
            
            
            
    def node_pf_colors(self, nodes, cell_stats):
        
            
        # zas, нет ножек - нет мультиков
        # node_pf_estimations nemusejí bejt
        try:
            # zkusmě pro jednoduchost 
            # čírou RGB hračku
            npf = nodes.node_pf_estimations
            cm = pg.colormap.ColorMap((0,1), [(0, 255, 0, 255), (255, 0, 0, 255)])
            return cm.mapToFloat(npf)
            
        except BaseException as e:
            msg = "node_pf_coloring has problems "
            error_msg = self.__class__.__name__ + ": " + msg + repr(e)
            print(error_msg)
            # simple coloring
            event = cell_stats['event']
            return self.get_color(event)
        

        
        
    
    def simple_coloring(self, nodes=None, cell_stats=None, *args, **kwargs):
        """
        if nodes and cell_stats  provided we will add them to self.cells
        otherwise function redraw items in self.cells
        """
        symbol_size = self.param.getValues()['node (pixel) size'][0]
        
        if nodes is None:
            for nodes, plot_item, cell_stats in self.cells:
                event = cell_stats['event']
                color = self.get_color(event)
                #symbolSize = np.sqrt(nodes.w / min(nodes.w)) # not bad
                plot_item.setData(color=color, size=symbol_size)
        
        # máme nodes, tj. jedeme poprvé        
        else:
            # draw tečičky
            #
            pos = getattr(nodes, self.sb_item.space)
            
            event = cell_stats['event']
            color = self.get_color(event)
            #symbolSize = np.sqrt(nodes.w / min(nodes.w)) # not bad
            plot_item = gl.GLScatterPlotItem(pos=pos, size=symbol_size, color=color)
            self.sb_item.central_widget.addItem(plot_item)
            
            # uložíme data
            self.cells.append((nodes, plot_item, cell_stats))
            
            # keep the GUI responsive :)
            self.sb_item.app.processEvents()
    
        
            
    def cell_probability_coloring(self, nodes=None, cell_stats=None, *args, **kwargs):
        """
        if nodes and cell_stats  provided we will add them to self.cells
        otherwise function redraw items in self.cells
        """
        symbol_size = self.param.getValues()['node (pixel) size'][0]
        
        if nodes is None:
            for nodes, plot_item, cell_stats in self.cells:
                event = cell_stats['event']
                cell_probability = cell_stats['cell_probability']
                if self.p_cell_max[event] < cell_probability:
                    self.p_cell_max[event] = cell_probability
                
            # přebarvíme tečičky podle obsahu pravděpodobnosti
            for nodes, plot_item, cell_stats in self.cells:
                # draw 
                pos = getattr(nodes, self.sb_item.space)
                #x, y = (*getattr(nodes, self.sb_item.space).T,)
                
                event = cell_stats['event']
                cell_probability = cell_stats['cell_probability']
                # bez modrého - maximální intenzita
                blue_intensity = 1 - cell_probability / self.p_cell_max[event]
                color = self.get_color(event, blue_intensity)
                #symbolSize = np.sqrt(nodes.w / min(nodes.w)) # not bad
                plot_item.setData(color=color, size=symbol_size)
        
        # máme nodes, tj. jedeme poprvé        
        else:
            event = cell_stats['event']
            cell_probability = cell_stats['cell_probability']
            # zkontrolujeme probability
            if self.p_cell_max[event] < cell_probability:
                self.p_cell_max[event] = cell_probability
                # a hned všecko dotyčné přebarvíme podle obsahu pravděpodobnosti
                for cell in self.cells:
                    if event == cell[2]['event']:
                        # bez modrého - maximální intenzita
                        # zde cell_probability se rovná self.p_cell_max[event]
                        # ale bacha! Nesplet se jinde!
                        blue_intensity = 1 - cell[2]['cell_probability'] / cell_probability
                        color = self.get_color(event, blue_intensity)
                        # bacha, potřebuji prvek vložit zpätky
                        cell[1].setData(color=color)
            
            # bez modrého - maximální intenzita
            blue_intensity = 1 - cell_probability / self.p_cell_max[event]
            color = self.get_color(event, blue_intensity)
            
            
            # draw tečičky
            #
            pos = getattr(nodes, self.sb_item.space)
            plot_item = gl.GLScatterPlotItem(pos=pos, size=symbol_size, color=color)
            self.sb_item.central_widget.addItem(plot_item)
            
            # uložíme data
            self.cells.append((nodes, plot_item, cell_stats))
            
            # keep the GUI responsive :)
            self.sb_item.app.processEvents()


    def get_color(self, event, blue_intensity=None):
        """
        get color for 'simple_coloring' or 'cell_probability_coloring'
        """
        # generally
        if event == 'success':
            color = [167, 255, 181]# xkcd:light seafoam green #a7ffb5
        elif event == 'failure':
            color = [253, 193, 197] # xkcd: pale rose (#fdc1c5)
        # já vím, že Voronoi nemá 'mix', эн но юа
        else:# 'mix'
            color = [255, 243, 154] # let's try xkcd: dark cream (#fff39a)
            
        if blue_intensity is not None:
            # let's play with blue a little bit
            # chcu mít korektní výstup
            # aspoň v něčem chcu být jistý
            if blue_intensity > 1:
                color[2] = 255
            elif blue_intensity < 0:
                color[2] = 0
            else:
                # pyqtgraph žere barvy i s čarkou
                # ale my je mu davat nebudeme
                color[2] = int(blue_intensity*255) 
            
        return (*(np.array(color)/255), 1)


"""
===========
♥ Чыры-пыры 
č Jiné
E Miscellaneous
===============
+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
"""

        
        
class CandidatesWidget(qt_plot.CandidatesWidget):


    def run_stm(self):
        #č indikace
        #self.setDisabled(True)
        with pg.BusyCursor():
            
            color_map = self.gradient.colorMap()
            
            try:#č může se tu stat cokoliv
                #č načíst sloupce prostě z libovolného vzorku
                cb = self.sb_item.sample_box.candidates_index.values().__iter__().__next__()
                #č je třeba mít seznam aktualní
                self.attr.setItems(list(cb.df.columns))
                
                #č neplest s self.attr!
                attr = self.attr.currentText()
                
                #č kruci, nejdřív je třeba najít maxvalue, minvalue je implicitně nula
                maxvalue = -np.inf
                minvalue = np.inf
                for id, cb in self.sb_item.sample_box.candidates_index.items():
                    array = getattr(cb, attr)
                    maxcb = np.nanmax(array)
                    mincb = np.nanmin(array)
                    if maxcb > maxvalue:
                        maxvalue = maxcb
                    if mincb < minvalue:
                        minvalue = mincb
                
                
                
                #č a teď jdeme!
                for id, cb in self.sb_item.sample_box.candidates_index.items():
                    array = getattr(cb, attr)
                    if np.isnan(array).any():
                        msg = "%s candidates has nans in %s attribute"%(id, attr)
                        error_msg = self.__class__.__name__ + ": " + msg
                        print(error_msg)
                    mask = np.isfinite(array)
                    values = array[mask]
                    norm_values = (values - minvalue) / (maxvalue - minvalue)
                    
                    pos = getattr(cb, self.sb_item.space)[mask]
                    color = color_map.mapToFloat(norm_values)
                    pen = gl.GLScatterPlotItem(pos=pos, size=self.sb_item.px_size, color=color)
                    self.sb_item.central_widget.addItem(pen)
                    self.pens.append(pen)
                
            except BaseException as e:
                msg = ""
                error_msg = self.__class__.__name__ + ": " + msg + repr(e)
                print(error_msg)
                self.error.emit(error_msg)
            
            
            
        # indikace
        #self.setEnabled(True)
        

    
        

        


Mode Type Size Ref File
100644 blob 17474 79f2258efa33670c57a70e41eff5d117bf69dcb7 IS_stat.py
100644 blob 6 0916b75b752887809bac2330f3de246c42c245cd __init__.py
100644 blob 63303 4a921c50676e2e8d1fd12a004d8dad7d2bcd1d38 blackbox.py
100644 blob 11078 cd8e7e5ab6fa5edde14878d8fa90bc65012ccc0d candybox.py
100644 blob 31020 1a4e4f3ccbfee2218309ffb33d1cf46bc249d4b9 estimation.py
100644 blob 18660 c1ef0f2dd7e992953ebd0b295de5c5c6721215d6 f_models.py
100644 blob 31025 70bab60405bfe783a2f7a9f2c41b7c1629d3d474 g_models.py
100644 blob 52315 cbecebc34eb21cc476b640b6c95ca9ba5ce745e1 gl_plot.py
100644 blob 2718 5d721d117448dbb96c554ea8f0e4651ffe9ac457 gp_plot.py
100644 blob 29393 96162a5d181b8307507ba2f44bafe984aa939163 lukiskon.py
100644 blob 10489 1f6dd06a036fdc4ba6a7e6d61ac0b84e8ad3a4c1 mplot.py
100644 blob 993 593f0411936cebc2b6ac12775f1b7c62e53e717c plot.py
100644 blob 2807 1feb1d43e90e027f35bbd0a6730ab18501cef63a plotly_plot.py
100644 blob 61709 379db4a6236832e448f387adb9063ab547e53e1f qt_plot.py
100644 blob 57475 645d9c4a2339f5a3997f1f16e61b58d0defba6cc qt_plot_old.py
100644 blob 6304 7fc6ac75e415df43af5b7aa9d6d1848aa5d0963d reader.py
100644 blob 4284 a0e0b4e593204ff6254f23a67652804db07800a6 samplebox.py
100644 blob 5553 bac994ae58f1df80c7f8b3f33955af5402f5a4f3 sball.py
100644 blob 2800 109ef02b48b0e735d9081be69bffd0471675bc98 simplex.py
100644 blob 21623 281aef80556b8d22842b8659f6f0b7dab0ad71af 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