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".
List of commits:
Subject Hash Author Date (UTC)
simplex: preliminary commit of an _Sense's auxiliary class - brand new SeparationAxis c0bbb714400df48abf7a398e018c5626803bf024 I am 2023-03-04 10:17:12
simplex._Sense: one more optimization a51f4ba65ad4a55ecc74fd82a6e011a4eeed8a5d I am 2023-02-26 15:21:23
qt_gui.qt_plot: implement sensitivity-related Arrows class b5365ae3ca3ab0ed8018e60d4e8711c5c30e4d96 I am 2023-02-26 14:31:10
simplex._Sense: rename "sensibility" to sensitivity a9d56d5386730263eb27e6e31f07b99bab0620fb I am 2023-02-26 14:29:54
g_models: add quadratic 013b4ddc108b94061eaebc1a10d18427f10a34d4 I am 2023-02-24 08:04:12
simplex._Sense: one more performance trick 59b423cca53b9975da67d876110927f233506de8 I am 2023-02-24 08:03:23
simplex: implement separability-based sensibility analysis (new brand _Sense class) 9c5d58f2301893ceaec1b0e90bff76035cfa15b2 I am 2023-02-23 18:49:11
dicebox.circumtri.CirQTri: switch to GaussCubatureIntegration 5b52dd25cb7a997af4953230116deb9efc577d56 I am 2023-02-11 04:32:48
simplex: implement GaussCubatureIntegration in the most memory-friendly way 689d253ae7e2a22242258fd5bef0a069caf7cf75 I am 2023-02-11 04:31:11
convex_hull.QHullCubature: implement memory-friendly outside (chi) integration ad8210a04b1e0903de7435cad16b1304707d0e6e I am 2023-02-09 22:22:05
qt_gui.qt_plot: require box recommendation to show next point 6f726047f7f08e884359020eaa1eac6f6cc125d2 I am 2023-02-09 11:51:44
dicebox.circumtri.CirQTri.get_circum_node: limit circumcenter by explore radius, not by just R 136ec73bb06da16c1f2bce64b3c349be4c8ba975 I am 2023-02-09 03:09:51
dicebox.circumtri: implement new brand CirQTri box 5879b8ad6317c65aa9481b59f76b6159f19e247a I am 2023-02-09 01:29:10
simplex.FullCubatureIntegration: store simplex probabilities in sorted dicts c0da90731ff3ede47d9b4eec0ad9b28a29027167 I am 2023-02-09 01:23:14
dicebox.circumtri: exploratory: even better idea 811ab11cd7172ff4a3807992f4928be2e8068ec0 I am 2023-02-08 15:31:23
dicebox.circumtri: exploratory, new idea 526d3f6185887ff48a341b0705d74dde4d15ca87 I am 2023-02-08 03:03:41
dicebox.circumtri: exploratory 806063d2e223c812280dc4845153450dd47faed3 I am 2023-02-06 17:15:15
dicebox.circumtri: exploratory efed2589f642d502f30e80f0e9b45dfeecd1c7c7 I am 2023-02-06 13:40:24
dicebox.circumtri: exploratory - again 34d3f4e47420e1c1e26b09570fb44d3270194751 I am 2023-02-06 12:50:45
qt_gui.qt_dicebox: change default q of circumtri classes 9fd5855e5d7cacf80d27fb383dd18a92d60e138b I am 2023-02-06 12:30:27
Commit c0bbb714400df48abf7a398e018c5626803bf024 - simplex: preliminary commit of an _Sense's auxiliary class - brand new SeparationAxis
Author: I am
Author date (UTC): 2023-03-04 10:17
Committer name: I am
Committer date (UTC): 2023-03-04 10:17
Parent(s): a51f4ba65ad4a55ecc74fd82a6e011a4eeed8a5d
Signer:
Signing key:
Signing status: N
Tree: bd620e14e70ed0e9d5346705a2ef0deee4606393
File Lines added Lines deleted
wellmet/simplex.py 239 0
File wellmet/simplex.py changed (mode: 100644) (index af2c23e..4d6d335)
... ... class _Sense:
542 542
543 543
544 544
545 #č tohle vůbec není žádná samostatná třída
546 #č Její jedina instance taky vůbec není samostatná.
547 #č simplex_scope je extrémně vázano na _Sense.process_simplex()
548 class SeparationAxis:
545 549
550 def __init__(self, sx, simplex_scope_reference):
551 #č simplex_scope potřebujeme pouze v jednom místě
552 self.simplex_scope = simplex_scope_reference
553 self.point_scope = set()
554 self.sx = sx # for sx.get_simplex_normal()
555
556 self.ED = points = sx.tri.points
557 self.failsi = failsi = sx.sample_box.failsi
558 self.simplices = sx.tri.simplices
559
560 nsim, nvar = points.shape
561
562 self.c_zeros = np.zeros(nvar + 1) # prices in terms of LP
563 self.b_ub = -np.ones(nsim)
564
565 self.A_ub = np.empty((nsim, nvar + 1))
566 self.A_ub[:, :-1] = points
567 self.A_ub[:, -1] = self.b_ub
568 self.A_ub *= (failsi * 2 - 1).reshape(-1, 1)
569
570
571
572 # method to start new task
573 def reset(self, simplex_id, finalized_nornal):
574 self.finalized_nornal = finalized_nornal
575
576 #č I když se mi nechce semka tahnout get_simplex_normal(),
577 #č stejně potřebuji indices pro point_scope
578 if finalized_nornal is None:
579 self.point_scope.clear()
580 indices = self.simplices[simplex_id]
581 self.point_scope.update(indices)
582 #č hloubka je nulová, žádné přípravné otačky konat nebudou.
583 #č tahle normála půjde rovnou na dějiště
584 #č poslední složka vektoru je buď nvar+1 dimenze
585 #č v případě get_simplex_normal()
586 #č nebo offset aka bias
587 #
588 # the last item of vector is either nvar+1 coordinate
589 # in case of get_simplex_normal()
590 # or offset aka bias
591 # in case of separation axis
592 # i.e. we should ignore it anyway
593 self._set_normal(self.sx.get_simplex_normal(indices)[:-1])
594
595
596
597
598
599 def _set_normal(self, normal):
600 #č podtržitko v názvu je k tomu, abych zdůraznil,
601 #č že point_scope musí být před volaním metody
602 #č jíž nastaven, musí být platný, musí být v pořádku
603 self.finalized_nornal = None
604 self.normal = normal
605
606
607 ED = self.ED
608 failsi = self.failsi
609
610 # set supports
611 max_red = -np.inf
612 min_green = np.inf
613 for point in self.point_scope:
614 failure = failsi[point]
615 value = np.inner(normal, ED[point])
616 if failure and (value > max_red):
617 max_red = value
618 self.red_supp = point
619 elif not failure and (value < min_green):
620 min_green = value
621 self.green_supp = point
622
623 assert max_red < min_green
624
625 self.max_red = max_red
626 self.min_green = min_green
627
628
629
630 def is_separable(self, simplex_id):
631 if self.finalized_nornal is None:
632 return self._is_non_finalized_separable(simplex_id)
633
634 if simplex_id in self.finalized_nornal.separated_simplices:
635 return True
636
637 #č předpokladáme, že finální vektor je omezen nadbytečnými simplexy
638 #č a tedy i nadbytečným vzálenými vzorky.
639 #č Že když pustíme LP my, tak kvůli menší bázi
640 #č zvládneme odseparovat i další simplexy.
641 #č To je poznámka k tomu,
642 #č že ze žádného porovnání množin simplexů a vzorků
643 #č nejsme shopní udělat žádný závěr.
644 #č Vyjímečný případ, že by množiny byly si rovny neuvažujeme.
645 #if self.finalized_nornal.separated_simplices == self.scope.keys():
646 # return False
647
648
649 # non-separable, then
650 #оӵ LP лэзьыны кулэ
651 #č poprvé s tomto runu potřebujeme tečky
652 #č odsuď jsou jenom dvě cesty:
653 #č ono se buď objeví být separabilní,
654 #č tehdy budeme pokračovat s nefinalizovaným vektorem.
655 #č nebo separabilní nebude, tehdy končíme, jdeme na další run.
656 simplices = self.simplices
657
658 point_scope = self.point_scope
659 #č vytvořme platný skoup
660 point_scope.clear()
661 #č jediné místo, kde simplex_scope potřebujeme!
662 for key in self.simplex_scope.keys():
663 point_scope.update(simplices[key])
664
665 #č máme finitně-finální vektor.
666 #č pokud teď čerstvě přídaný simplex separovat nepůjde,
667 #č tak volající kód timhle skončí, ten finitní vektor vrátí.
668 #č Když půjde, tak budeme pokračovat s nově nalezenou normálou.
669 #č To je k tomu, že nemá cenu brat v uvahu první scenář.
670 #č Tam to hned skončí a je to jedno, co tu budeme mít za skoupy.
671 #č A nejen to. Kdyby někdo neohledúplně volal is_separable(),
672 #č tak v případě False my nesaháme na finalizovaný vektor
673 #č a tak v příštím volaní point_scope se očístí a naplní se znovu.
674 #č Mimochodem, díky takovému odvažnému apdejtu je nám jedno,
675 #č jak přesně volající kód zachází s simplex_scope,
676 #č zda přídává simplex před, po, zárověň...
677 point_scope.update(simplices[simplex_id])
678
679 result = self.get_separation_axis(list(point_scope))
680 if not result.success: #č konec, zde končíme
681 if result.status != 2: # Problem appears to be infeasible.
682 print("Sense: linprog ended with status %s" % result.status)
683 return False #č na finalized_nornal nesáháme
684
685 # sucсess!
686 #č musíme se přípravit k dalšímu životu
687 #č _set_normal() tu totalní normalu vynulue
688 #self.finalized_nornal = None
689
690 #č hned odřízneme poslední složku -
691 #č bude tam posunutí b aka offcet aka bias
692 #č point_scope jíž máme v pořádku.
693 self._set_normal(result.x[:-1])
694
695 return True
696
697
698
699 #č Pokud zkoušíme separabilitu sousedícího simplexu,
700 #č tak narazíme na maximálně jeden nový vzorek.
701 #č Zbytek nového simplexu tvoří společná stěna-hyperrovina.
702 def _is_non_finalized_separable(self, simplex_id):
703
704 #č tady je to opakový opak
705 #č Pokud vzorek nejde separovat, tak se nesmí dostat do point_scope
706 #č V každém případě se bude pokračovat až do posledního simplexu
707
708 #č metoda buď přída vzorek do skoupu a vratí None,
709 #č nebo nepřídá a ho vrátí. Striktně jeden.
710 non_separable_point = self._is_neighbor_separable(simplex_id)
711 if non_separable_point is None:
712 return True
713
714 #č point_scope už máme
715 point_list = list(self.point_scope)
716 point_list.add(non_separable_point)
717
718 result = self.get_separation_axis(point_list)
719 if not result.success: #č konec, zde končíme
720 if result.status != 2: # Problem appears to be infeasible.
721 print("Sense: linprog ended with status %s" % result.status)
722 return False
723
724 # sucсess!
725 #č dodáme ten diskutovaný pochybný, jíž ospravedlněný vzorek
726 self.point_scope.add(non_separable_point)
727
728 #č hned odřízneme poslední složku -
729 #č bude tam posunutí b aka offcet aka bias
730 self._set_normal(result.x[:-1])
731
732 return True
733
734
735
736 #č Pokud zkoušíme separabilitu sousedícího simplexu,
737 #č tak narazíme na maximálně jeden nový vzorek.
738 #č Zbytek nového simplexu tvoří společná stěna-hyperrovina.
739 #č zkusme takhle. V případě, že tamten nový vzorek jde odseparovat,
740 #č tak vratíme None, jinak číslo toho vzorku pro další kód
741 def _is_neighbor_separable(self, simplex_id):
742 #č z logiky volajícího kódu nema cenu kontrolovat simplex_scope
743 #č ono se ptá na separabilitu pouze toho, co jíž není ve skoupu
744 #if simplex_id in self.simplex_scope:
745 # return True
746
747 for point in self.simplices[simplex_id]:
748 if point not in self.point_scope:
749 value = np.inner(self.normal, self.ED[point])
750 if self.failsi[point]:
751 if value > self.min_green:
752 return point
753 elif value > self.max_red:
754 self.max_red = value
755 self.red_supp = point
756 else:
757 if value < self.max_red:
758 return point
759 elif value < self.min_green:
760 self.min_green = value
761 self.green_supp = point
762
763 #č v této fázi,
764 #č když jsme doposud neskončíli,
765 #č už víme, že vrchol je separabilní
766 self.point_scope.add(point)
767
768 #č takovej vzorek, který není ve skoupu
769 #č u sousedicího simplexu může být maximálně jeden
770 return None
771
772
773 #č přes liblinear to taky jde
774 #from liblinear import liblinearutil
775 #pr = liblinearutil.problem(box.failsi, box.G)
776 #param = liblinearutil.parameter('-s 2 -c 10000 -B 1')
777 #m = liblinearutil.train(pr, param)
778 #m.get_decfun()
779 def get_separation_axis(self, indices):
780 A_ub = self.A_ub[indices]
781 b_ub = self.b_ub[:len(A_ub)]
782 return optimize.linprog(self.c_zeros, A_ub=A_ub, b_ub=b_ub,
783 options={'presolve':False}, method='highs-ds',
784 bounds=(None, None))
546 785
547 786
548 787
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