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._Triangulation: implement averaged gradients() f6608597cf490dd6cf1a6af3b7d8e6d0a72732d9 I am 2023-01-15 09:37:14
add line and two lines problems for sensitivity analysis 0c620f7960399ed8b34e02d663772d0eb34e66f1 I am 2023-01-14 08:30:56
simplex.FullCubatureIntegration.get_sensitivities: do not normalize global gradient c648c484a5599a00236ab70c36cdf407c24f5306 I am 2023-01-14 04:51:16
simplex: implement sensitivities 65056319c51fd93deb76ce14e52c08869dfc04a2 I am 2023-01-13 11:27:46
mplot.maxes: prepare GRaph plot 68774535915af3deaf6104a0bbe9b3c4d24c0cec I am 2023-01-12 02:43:20
mplot.mgraph.tri_estimation_plot: use planar vertex estimator 1d300fa56248be9cdab5a4e5b0d4372f68fd5255 I am 2023-01-11 15:48:28
mplot.mart: implement lsf_boundary() contouring 8642d2f0c043f4e71af1cc45b9e25da326afd3ca I am 2023-01-11 05:20:00
mplot.mart.setup_labels: make labels more adaptive 60cd0b34a67b38df1c730dddc9e0771139818f33 I am 2023-01-11 04:35:47
simplex._SamplingTriangulation: self reference fix 79204ff15cba7c5ca24bb5b48e2073058acf00ae I am 2023-01-11 03:21:34
testcases.testcases_2D: add pf for rastrigin 4f3b35fe18496fcb9f59860174c876cb7a82749a I am 2023-01-10 23:39:45
dicebox.circumtri: fix initialization on ED loading. Use simple potential for space filling c02656e588a164685c4a750ce865c4b50f716563 I am 2023-01-10 23:37:30
simplex.FullCubatureIntegration.get_failure_moments: negative probabilities fix 27dd6bb1e5bb0eb1e28eef88e28913d1eca461da I am 2023-01-10 20:16:09
add more nD problems from literature 2380cdbdcf17dd5e0af19506c6f7bb1e838de03b I am 2023-01-10 06:13:56
f_models: implement Norm class 698c3216b1c8a0e2d57bfb33340a7ea8758d4e63 I am 2023-01-09 23:13:48
simplex.FullCabatureIntegration: implement get_failure_moments() 0742ee3d3932d7c3c0b6bbe6e327c92ced04aedd I am 2023-01-09 19:53:55
simplex._Triangulation: add get_failure_moments() method 3e362970b67a50643c1e486729efbf554cacf03f I am 2023-01-08 03:26:02
simplex: use advanced indexing instead of isin function 568c8623ad4cfaa343df05ae880d099556037cf0 I am 2023-01-06 13:41:00
dicebox.circumtri: fill in global simplex index only in case of full integration 065668f1810396f55be0860416f1c6ad507a9573 I am 2023-01-04 21:44:57
dicebox.circumtri: implement "holyday", breaks the greatest safe simplex 14f067d43cde14f990e7a4136b61c548f3acde92 I am 2023-01-04 16:28:51
simplex: implement FullCubatureIntegration 072927115a8a641f89b39265085eb6cac27aa492 I am 2023-01-04 16:27:39
Commit f6608597cf490dd6cf1a6af3b7d8e6d0a72732d9 - simplex._Triangulation: implement averaged gradients()
Author: I am
Author date (UTC): 2023-01-15 09:37
Committer name: I am
Committer date (UTC): 2023-01-15 09:37
Parent(s): 0c620f7960399ed8b34e02d663772d0eb34e66f1
Signer:
Signing key:
Signing status: N
Tree: 4fabf9fc93b559bd87655c2aba9d12ca9bd8a3bf
File Lines added Lines deleted
wellmet/simplex.py 55 3
File wellmet/simplex.py changed (mode: 100644) (index 979004b..471a2db)
... ... class _Triangulation:
399 399 return np.logical_xor(has_failure, all_failure) return np.logical_xor(has_failure, all_failure)
400 400
401 401
402 def get_mixed_normals(sx):
403 simplices = sx.tri.simplices[sx.is_mixed(sx.tri.simplices)]
402 def get_averaged_mixed_normals(sx):
403 simplices = sx.tri.simplices
404
405 nvar = sx.sample_box.nvar
406 failsi = sx.sample_box.failsi
407 mask = np.empty(sx.tri.nsimplex + 1, dtype=bool)
408 mask[-1] = False
409
410 in_failure = failsi[simplices]
411 has_failure = in_failure.any(axis=1)
412 all_failure = in_failure.all(axis=1)
413 np.logical_xor(has_failure, all_failure, out=mask[:-1])
414
415
416 mixed_simplices = simplices[mask[:-1]]
417
418
419
420 X = np.empty((sx.sample_box.nsim, nvar + 1))
421 X[:, :-1] = getattr(sx.sample_box, sx.tri_space)
422 X[:, -1] = failsi
423
424 # X = np.empty((nvar + 1, sx.sample_box.nsim))
425 # X[:-1] = getattr(sx.sample_box, sx.tri_space).T
426 # X[-1] = failsi
427
428
429 vectors = X[mixed_simplices[:, 1:]] - X[mixed_simplices[:, :1]]
430 basises, __ = np.linalg.qr(vectors.transpose(0, 2, 1), mode='complete')
431
432 #č o jeden simplex víc aby index -1 odkazoval na nulovou normálu
433 normals = np.zeros((sx.tri.nsimplex + 1, nvar + 1))
434 normals[mask] = basises[:, :, -1]
435
436 #č ujistit se, že normály směrovany nahoru,
437 #č tj. směrem k bezpečným vzorkám
438 normals *= (1 - 2 * (normals[:,-1:] < 0))
439
440 neighbors = sx.tri.neighbors[mask[:-1]]
441
442 mixed_normals = normals[mask]
443 # n_mixed_simplices x nvar+1 (neighbors) x nvar+1 (coordinates)
444 mixed_normals += np.sum(normals[neighbors], axis=1)
445
446 return mixed_normals, mixed_simplices
447
448
449 def get_pure_mixed_normals(sx):
450 simplices = sx.tri.simplices
451
452 mask = sx.is_mixed(simplices)
453 simplices = simplices[mask]
454
404 455 nvar = sx.sample_box.nvar nvar = sx.sample_box.nvar
405 456 failsi = sx.sample_box.failsi failsi = sx.sample_box.failsi
406 457
 
... ... class _Triangulation:
426 477
427 478
428 479 def get_gradients(sx): def get_gradients(sx):
429 normals, simplices = sx.get_mixed_normals()
480 normals, simplices = sx.get_averaged_mixed_normals()
430 481
431 482 gradients = normals[:, :-1] gradients = normals[:, :-1]
432 483
 
... ... class FullCubatureIntegration(_Triangulation):
1399 1450 #č odhady jsou ve slovníku, posíláme jen to, co tam není #č odhady jsou ve slovníku, posíláme jen to, co tam není
1400 1451 sx.on_failure_added(simplex, indices, vertices_model, nodes, vol) sx.on_failure_added(simplex, indices, vertices_model, nodes, vol)
1401 1452
1453
1402 1454 def get_sensitivities(sx): def get_sensitivities(sx):
1403 1455 gradients, simplices = sx.get_gradients() gradients, simplices = sx.get_gradients()
1404 1456
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