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)
introduce shell module, clean version of ghull b13e9b2b15d109770e0c3fbdd5bb7b2a3b5bf5fc I am 2022-12-28 15:45:27
qt_gui.qt_plot.CandidatesWidget: regression fix a75874b5a53e7b420672b582f3ff415b17a32e9b I am 2022-12-23 12:19:51
qt_gui.qt_pairwise.ContactWidget: add mask function 6f88fec7db92e5cd21826a06f42a609a3e03cb6b I am 2022-12-21 15:20:01
qt_gui.qt_pairwise.ContactWidget: add option for mixed pairs only adjacency search. 2f87e8a18908de607b8abc0779448ca3095c134c I am 2022-12-21 14:30:23
qt_gui.qt_dicebox: reflect changes in CircumTri class 9821d389240a383f0b13b7dd92cd501e1c40c537 I am 2022-12-15 09:36:16
dicebox.circumtri: employ _Exploration class. Keep naive CircumTri class in renamed module b8605fb24d55b3296ceac34afba541a5b7ee3070 I am 2022-12-15 09:33:59
convex_hull.QHull.get_exploration_vector: self -> hull fix dae47ef37b5871a98d2ea274cc35aeb2a4eebe1c I am 2022-12-15 09:30:31
qt_gui/qt_dicebox: add DumpExplorationWidget b687895f06c780418a12eea59bd052ad7cd348f3 I am 2022-12-15 01:31:14
dicebox: add helper _exploration module 836ae6cd12f1a04f1e2fc3423e46777d16ebad7d I am 2022-12-15 01:30:20
convex_hull: add get_exploration_vector() method 795259d91b64c732e60a83a1d9e3b89c2141b003 I am 2022-12-15 01:28:12
qt_gui.qt_gui: dot not create graph widgets if estimations are not present f070ed140a135a980e4436ed3779455c0e1bf399 I am 2022-12-14 02:25:45
qt_gui.gl_plot.CandidatesWidget: apply changes for GLplot too 03550040b60c67d247e26425e35fac565dd88886 Aleksei Gerasimov 2022-12-14 01:56:01
qt_gui.qt_plot.CandidatesWidget: add CandyNodes support 2357cdd92170b1f660b15c92aec695d3f7c108b3 I am 2022-12-14 01:08:30
qt_gui.qt_graph_widgets.EstimationGraph: reimplement piece of zerosafe logic 342a1b7bfddb7ba5695b3f0e017c7f75811db1f0 I am 2022-12-13 22:15:54
qt_gui.qt_graph_widgets: update x range when box runned c52052a9c3f6926d03e8c39f7642aa0a8864fdc3 I am 2022-12-13 10:23:45
qt_gui.qt_graph_widgets.EstimationGraph: add outside and mixed lines df6c06e4287492855d613b1ae46edfca6a4757d3 I am 2022-12-12 06:13:55
qt_gui: replace old copypaste by cleaned up EstimationGraph 515d0d8a892984f8a7667da874972daff3318d74 I am 2022-12-12 05:36:43
qt_gui: prepare new BoxEstimationData class; rework ErrorGraph 37387100715f980771a8a37f99eec3e723759cd3 I am 2022-12-11 23:17:06
dicebox.circumtri: replace ndim by nvar to make TriEstimation pandas-friendly 0824222ae4c9c1f3fede971481816ed9ded9707a I am 2022-12-11 23:14:05
step back: introduce box_estimations, keep "estimations" for user ones dde983d57423e28083980cbae05f0a81fc8311e6 I am 2022-12-10 03:54:09
Commit b13e9b2b15d109770e0c3fbdd5bb7b2a3b5bf5fc - introduce shell module, clean version of ghull
Author: I am
Author date (UTC): 2022-12-28 15:45
Committer name: I am
Committer date (UTC): 2022-12-28 15:45
Parent(s): a75874b5a53e7b420672b582f3ff415b17a32e9b
Signer:
Signing key:
Signing status: N
Tree: c129c0167def07ec01226a0612640a1f90434a26
File Lines added Lines deleted
wellmet/shell.py 238 0
File wellmet/shell.py added (mode: 100644) (index 0000000..d4c45fd)
1 #!/usr/bin/env python
2 # coding: utf-8
3
4 import numpy as np
5 from . import sball
6 from scipy import stats
7 from .candynodes import CandyNodes
8 from .IS_stat import get_1DS_sample # for Shell_1DS
9
10 from collections import namedtuple
11
12 try: # try to outthink OS's memory management
13 import os
14 mem_bytes = os.sysconf('SC_PAGE_SIZE') * os.sysconf('SC_PHYS_PAGES')
15 except BaseException as e:
16 mem_GB = 16
17 mem_bytes = mem_GB * 1024**3 # hello, Windows!
18 #print("ghull failed to get amount of RAM installed. %s GB is supposed."% mem_GB, repr(e))
19 #print("BTW, you are using Windows, aren't you?")
20
21
22
23 ShellStats = namedtuple('ShellStats', (
24 'nsim',
25 'nvar',
26 'nfacets',
27 'r', 'R',
28 'inner',
29 'shell',
30 'outer',
31 'FORM_outside',
32 'TwoFORM_outside',
33 'orth_outside'
34 ))
35
36 ShellEstimation = namedtuple('ShellEstimation',
37 ('shell_budget', 'shell_inside',
38 'shell_outside', 'inside', 'outside'))
39
40
41 #č pomocná třída pro integrování
42 #č nejsem tedy jist, jestli je to nejlepší napad -
43 #č dělit Ghull na podtřídy, delegovat funkcionalitu
44 #č ale jinak Ghull se stavá hodně překomplikovaným.
45 #č nic lepšího mně nenapadá, přemyšlel jsem dlouho.
46 class Shell_1DS:
47 """1DS stands for 1D sampling.
48 1DS is, actually, an importance sampling method
49 that, actually, does not calculate IS weights
50 as f_density / h_density, but
51 (arbitrary, nonuniformly) divides interval to subintervals.
52 By using CDF transformes subintervals to 0-1 measure line.
53 Then gets one node as middle point of every subinterval,
54 weights therefore are just interval widths itself.
55 No sampling imprecisions are introduced,
56 therefore no spring, no correction are needed.
57 """
58 def __init__(self, hull, shell):
59 self.shell = shell
60 self.hull = hull
61 self.nvar = hull.sample.nvar
62
63 #č objevilo se, že scipy.stats.chi je neunosně nepřesné
64 #č vzdává se jíž na poloměru 8
65 self.sball = sball.Sball(self.nvar)
66
67 self.integration_cutoff = np.inf
68
69
70 def integrate(self, nis, callback_all=None, callback_outside=None):
71 self.reset(nis)
72
73 if self.hull.nsimplex == 0:
74 bus = self.integration_cutoff
75 else:
76 bus = int(mem_bytes / self.hull.nsimplex / 8 / 10) + 1
77 while self.nsampled < nis:
78
79 seats = min(nis - self.nsampled, self.integration_cutoff, bus)
80 try:
81 self._process_part(seats, nis, callback_all, callback_outside)
82 except MemoryError as m:
83 print(self.__class__.__name__ +":", "memory error, %s sedaček" % seats, repr(m))
84 self.integration_cutoff = int(np.ceil(seats/2))
85
86 assert nis == self.nsampled
87
88 return self._get_result()
89
90
91 def reset(self, nis): # clear
92
93 self.nsampled = 0
94
95 #č poloměry bereme ze skořapky
96 #č za správné nastavení (stejně jako u MC)
97 #č zodpovidá uživatel třídy
98 #č třída vlastní odhad r nijak nevyuživá!
99 r = self.shell.r
100 R = self.shell.R
101 # let's predefine 1D sequence at very beginning
102 if r > np.sqrt(self.nvar - 1):
103 x_sub = np.geomspace(r, R, nis+1, endpoint=True)
104 else:
105 x_sub = np.linspace(r, R, nis+1, endpoint=True)
106
107
108 x, weights = get_1DS_sample(self.sball, x_sub)
109 self.x = x
110 self.weights = weights
111 self.mask = np.empty(nis, dtype=np.bool)
112
113
114
115 # bus analogy
116 def _process_part(self, seats, nis, callback_all=None, callback_outside=None):
117 # boarding
118 left = self.nsampled
119 right = self.nsampled + seats
120 rs = self.x[left:right]
121
122 # rand_dir: prepare ns random directions on a unit d-sphere
123 rand_dir = sball.get_random_directions(seats, self.nvar) #random directions
124 nodes_G = rand_dir*rs[:,None]
125 nodes = self.hull.sample.f_model.new_sample(nodes_G, space='G')
126
127 # who is back?
128 d = self.hull.get_hyperplane_distances(nodes)
129 mask = d > 0
130 assert len(mask) == seats
131
132 # the most important part
133 self.nsampled += len(mask)
134 self.mask[left:right] = mask
135
136 if callback_all is not None:
137 # -2 = 'inside' -1 = 'outside'
138 candy_nodes = CandyNodes(nodes, event_id=mask-2, is_outside=mask, d=d)
139 callback_all(candy_nodes)
140
141 if (callback_outside is not None) and np.any(mask):
142 callback_outside(nodes[mask], d[mask])
143
144
145
146 def _get_result(self):
147 #č mask related to hull.is_outside()
148 shell_pf = np.sum(self.weights[self.mask])
149 shell_ps = np.sum(self.weights[~self.mask])
150
151 return self.nsampled, shell_ps, shell_pf
152
153
154
155
156 class GaussianAnnulus:
157 def __init__(self, hull):
158 assert hull.space == 'G'
159
160 self.hull = hull
161 self.shell = sball.Shell(hull.sample.nvar)
162 self.outside_dist = sball.Shell(hull.sample.nvar)
163 self.sample = hull.sample
164
165 self.gint = Shell_1DS(self.hull, self.shell)
166
167
168 def boom(self, ns, use_MC=False):
169 if use_MC:
170 self.outside_dist.set_bounds(self.get_R())
171 nodes_G = self.outside_dist.rvs(ns)
172 else:
173 # rand_dir: prepare ns random directions on a unit d-sphere
174 rand_dir = sball.get_random_directions(ns, self.sample.nvar) #random directions
175
176 #č deme od vnější .get_R() kružnici směrem ven
177 r = self.get_R()
178
179 # maximum radius, where norm.pdf() wasn't zero
180 # -38.575500173381374935388521407730877399444580
181 # don't ask me what the magic python use to distinguish
182 # digits after double precision
183 max_R_ever = 37
184 if r < max_R_ever:
185 R = max_R_ever
186 else:
187 R = r + 10
188 r = np.linspace(self.get_R(), max_R_ever, ns, endpoint=True)
189 nodes_G = rand_dir*r[:,None]
190
191 nodes = self.sample.f_model.new_sample(nodes_G, space='G')
192 return nodes
193
194
195 def get_R(self):
196 sum_squared = np.sum(np.square(self.sample.G), axis=1)
197 #index = np.argmax(sum_squared)
198 return np.sqrt(np.nanmax(sum_squared))
199
200 def setup_shell(self):
201 shell = self.shell
202 r = self.hull.get_r()
203 R = self.get_R()
204
205 if r<0:
206 shell.set_bounds(0, R)
207 else:
208 shell.set_bounds(r, R)
209
210 return r, R
211
212 def get_shell_estimation(self):
213 r, R = self.setup_shell()
214
215 return ShellStats(
216 self.sample.nsim,
217 self.sample.nvar,
218 self.hull.nsimplex,
219 r, R,
220 self.shell.ps,
221 self.shell.p_shell,
222 self.shell.pf,
223 stats.norm.sf(r),
224 self.hull.get_2FORM_outside(),
225 self.hull.get_orth_outside()
226 )
227
228 def integrate(self, nis, callback_all=None, callback_outside=None):
229 self.setup_shell()
230
231 nsampled, shell_ps, shell_pf = self.gint.integrate(nis, callback_all, callback_outside)
232
233 inside = shell_ps + self.shell.ps
234 outside = shell_pf + self.shell.pf
235
236 return ShellEstimation(nsampled, shell_ps, shell_pf, inside, outside)
237
238
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