|
1 | 1 | from typing import List |
2 | 2 | import numpy as np |
3 | 3 |
|
4 | | -from dataclasses import dataclass |
| 4 | +from dataclasses import dataclass, field |
5 | 5 | from sectionproperties.pre.pre import Material |
6 | 6 |
|
7 | 7 |
|
@@ -36,6 +36,25 @@ class Tri6: |
36 | 36 | coords: np.ndarray |
37 | 37 | node_ids: List[int] |
38 | 38 | material: Material |
| 39 | + # _M and _x0 are used for global coord to local coord mapping |
| 40 | + _M: np.ndarray = field(init=False) |
| 41 | + _x0: np.ndarray = field(init=False) |
| 42 | + |
| 43 | + def __post_init__(self): |
| 44 | + # Create a mapping from global elm to local elm (unit triangle) |
| 45 | + # The result is used in the global to local mapping: (eta, xi) = _M(x_global-_x0), zeta = 1-eta-xi |
| 46 | + (p0, p1, self._x0) = self.coords[:, 0:3].transpose() |
| 47 | + # Shift the triangle so the x_3 vertex (global) is on the origin ((eta, xi, zeta) = (0,0,1)) |
| 48 | + # Aside: This is chosen to be consistent with the shape function definitions. |
| 49 | + # At (0,0,1), N3=1, N1=N2=N4=...=0 (see Theoretical Background documentation) |
| 50 | + # since (x, y) = (sum(N_i(eta,xi,zeta)*x_i), sum(N_i(eta,xi,zeta)*y_i) |
| 51 | + # then at (0,0,1): (x,y) = (x_3,y_3). ie: (x_3, y_3) => (0,0,1) |
| 52 | + r0 = p0 - self._x0 |
| 53 | + r1 = p1 - self._x0 |
| 54 | + # Asseble the equations to solve for the transformation for the unit triangle |
| 55 | + x = np.array([r0, r1]).transpose() |
| 56 | + b = np.array([[1, 0], [0, 1]]) |
| 57 | + self._M = np.linalg.solve(x, b) |
39 | 58 |
|
40 | 59 | def __repr__(self): |
41 | 60 | rep = f"el_id: {self.el_id}\ncoords: {self.coords}\nnode_ids: {self.node_ids}\nmaterial: {self.material}" |
@@ -536,6 +555,162 @@ def element_stress( |
536 | 555 | gps[:, 0], |
537 | 556 | ) |
538 | 557 |
|
| 558 | + def local_element_stress( |
| 559 | + self, |
| 560 | + p, |
| 561 | + N, |
| 562 | + Mxx, |
| 563 | + Myy, |
| 564 | + M11, |
| 565 | + M22, |
| 566 | + Mzz, |
| 567 | + Vx, |
| 568 | + Vy, |
| 569 | + ea, |
| 570 | + cx, |
| 571 | + cy, |
| 572 | + ixx, |
| 573 | + iyy, |
| 574 | + ixy, |
| 575 | + i11, |
| 576 | + i22, |
| 577 | + phi, |
| 578 | + j, |
| 579 | + nu, |
| 580 | + omega, |
| 581 | + psi_shear, |
| 582 | + phi_shear, |
| 583 | + Delta_s, |
| 584 | + ): |
| 585 | + """Calculates the stress at a point `p` within the element resulting from a specified loading. |
| 586 | +
|
| 587 | + :param p: Point (x,y) in the global coordinate system that is within the element. |
| 588 | + :type p: :class:`numpy.ndarray` |
| 589 | + :param float N: Axial force |
| 590 | + :param float Mxx: Bending moment about the centroidal xx-axis |
| 591 | + :param float Myy: Bending moment about the centroidal yy-axis |
| 592 | + :param float M11: Bending moment about the centroidal 11-axis |
| 593 | + :param float M22: Bending moment about the centroidal 22-axis |
| 594 | + :param float Mzz: Torsion moment about the centroidal zz-axis |
| 595 | + :param float Vx: Shear force acting in the x-direction |
| 596 | + :param float Vy: Shear force acting in the y-direction |
| 597 | + :param float ea: Modulus weighted area |
| 598 | + :param float cx: x position of the elastic centroid |
| 599 | + :param float cy: y position of the elastic centroid |
| 600 | + :param float ixx: Second moment of area about the centroidal x-axis |
| 601 | + :param float iyy: Second moment of area about the centroidal y-axis |
| 602 | + :param float ixy: Second moment of area about the centroidal xy-axis |
| 603 | + :param float i11: Second moment of area about the principal 11-axis |
| 604 | + :param float i22: Second moment of area about the principal 22-axis |
| 605 | + :param float phi: Principal bending axis angle |
| 606 | + :param float j: St. Venant torsion constant |
| 607 | + :param float nu: Effective Poisson's ratio for the cross-section |
| 608 | + :param omega: Values of the warping function at the element nodes |
| 609 | + :type omega: :class:`numpy.ndarray` |
| 610 | + :param psi_shear: Values of the psi shear function at the element nodes |
| 611 | + :type psi_shear: :class:`numpy.ndarray` |
| 612 | + :param phi_shear: Values of the phi shear function at the element nodes |
| 613 | + :type phi_shear: :class:`numpy.ndarray` |
| 614 | + :param float Delta_s: Cross-section shear factor |
| 615 | + :return: Tuple containing stress values at point `p` |
| 616 | + (:math:`\sigma_{zz,n}`, :math:`\sigma_{zz,mxx}`, |
| 617 | + :math:`\sigma_{zz,myy}`, :math:`\sigma_{zz,m11}`, |
| 618 | + :math:`\sigma_{zz,m22}`, :math:`\sigma_{zx,mzz}`, |
| 619 | + :math:`\sigma_{zy,mzz}`, :math:`\sigma_{zx,vx}`, |
| 620 | + :math:`\sigma_{zy,vx}`, :math:`\sigma_{zx,vy}`, |
| 621 | + :math:`\sigma_{zy,vy}`) |
| 622 | + :rtype: tuple(float, float, ...) |
| 623 | + """ |
| 624 | + |
| 625 | + # get the elements nodal stress |
| 626 | + ( |
| 627 | + sig_zz_n_el, |
| 628 | + sig_zz_mxx_el, |
| 629 | + sig_zz_myy_el, |
| 630 | + sig_zz_m11_el, |
| 631 | + sig_zz_m22_el, |
| 632 | + sig_zx_mzz_el, |
| 633 | + sig_zy_mzz_el, |
| 634 | + sig_zx_vx_el, |
| 635 | + sig_zy_vx_el, |
| 636 | + sig_zx_vy_el, |
| 637 | + sig_zy_vy_el, |
| 638 | + weights, |
| 639 | + ) = self.element_stress( |
| 640 | + N, |
| 641 | + Mxx, |
| 642 | + Myy, |
| 643 | + M11, |
| 644 | + M22, |
| 645 | + Mzz, |
| 646 | + Vx, |
| 647 | + Vy, |
| 648 | + ea, |
| 649 | + cx, |
| 650 | + cy, |
| 651 | + ixx, |
| 652 | + iyy, |
| 653 | + ixy, |
| 654 | + i11, |
| 655 | + i22, |
| 656 | + phi, |
| 657 | + j, |
| 658 | + nu, |
| 659 | + omega, |
| 660 | + psi_shear, |
| 661 | + phi_shear, |
| 662 | + Delta_s, |
| 663 | + ) |
| 664 | + |
| 665 | + # get the local coordinates of the point within the reference element |
| 666 | + p_local = self.local_coord(p) |
| 667 | + # get the value of the basis functions at p_local |
| 668 | + N = shape_function_only(p_local) |
| 669 | + |
| 670 | + # interpolate the nodal values to p_local and add the results |
| 671 | + ( |
| 672 | + sig_zz_n_p, |
| 673 | + sig_zz_mxx_p, |
| 674 | + sig_zz_myy_p, |
| 675 | + sig_zz_m11_p, |
| 676 | + sig_zz_m22_p, |
| 677 | + sig_zx_mzz_p, |
| 678 | + sig_zy_mzz_p, |
| 679 | + sig_zx_vx_p, |
| 680 | + sig_zy_vx_p, |
| 681 | + sig_zx_vy_p, |
| 682 | + sig_zy_vy_p, |
| 683 | + ) = np.dot( |
| 684 | + np.array( |
| 685 | + [ |
| 686 | + sig_zz_n_el, |
| 687 | + sig_zz_mxx_el, |
| 688 | + sig_zz_myy_el, |
| 689 | + sig_zz_m11_el, |
| 690 | + sig_zz_m22_el, |
| 691 | + sig_zx_mzz_el, |
| 692 | + sig_zy_mzz_el, |
| 693 | + sig_zx_vx_el, |
| 694 | + sig_zy_vx_el, |
| 695 | + sig_zx_vy_el, |
| 696 | + sig_zy_vy_el, |
| 697 | + ] |
| 698 | + ), |
| 699 | + N, |
| 700 | + ) |
| 701 | + sig_zz_m_p = sig_zz_mxx_p + sig_zz_myy_p + sig_zz_m11_p + sig_zz_m22_p |
| 702 | + sig_zx_v_p = sig_zx_vx_p + sig_zx_vy_p |
| 703 | + sig_zy_v_p = sig_zy_vx_p + sig_zy_vy_p |
| 704 | + sig_zy_p = sig_zy_mzz_p + sig_zy_v_p |
| 705 | + sig_zz_p = sig_zz_n_p + sig_zz_m_p |
| 706 | + sig_zx_p = sig_zx_mzz_p + sig_zx_v_p |
| 707 | + sig_zy_p = sig_zy_mzz_p + sig_zy_v_p |
| 708 | + return ( |
| 709 | + sig_zz_p, |
| 710 | + sig_zx_p, |
| 711 | + sig_zy_p, |
| 712 | + ) |
| 713 | + |
539 | 714 | def point_within_element(self, pt): |
540 | 715 | """Determines whether a point lies within the current element. |
541 | 716 |
|
@@ -571,6 +746,19 @@ def point_within_element(self, pt): |
571 | 746 | else: |
572 | 747 | return False |
573 | 748 |
|
| 749 | + def local_coord(self, p): |
| 750 | + """Map a point `p` = (x, y) in the global coordinate system onto a |
| 751 | + point (eta, xi, zeta) in the local coordinate system. |
| 752 | +
|
| 753 | + :param p: Global coordinate (x,y) |
| 754 | + :type p: :class:`numpy.ndarray` |
| 755 | + :return: Point in local coordinate (eta, xi, zeta) |
| 756 | + :rtype: :class:`numpy.ndarray` |
| 757 | + """ |
| 758 | + (eta, xi) = np.dot(self._M, p - self._x0) |
| 759 | + zeta = 1 - eta - xi |
| 760 | + return np.array([eta, xi, zeta]) |
| 761 | + |
574 | 762 |
|
575 | 763 | def gauss_points(n): |
576 | 764 | """Returns the Gaussian weights and locations for *n* point Gaussian integration of a quadratic |
@@ -675,6 +863,30 @@ def shape_function(coords, gauss_point): |
675 | 863 | return (N, B, j) |
676 | 864 |
|
677 | 865 |
|
| 866 | +def shape_function_only(p): |
| 867 | + """The values of the Tri6 shape function at a point `p`. |
| 868 | +
|
| 869 | + :param p: Point (eta, xi, zeta) in the local coordinate system. |
| 870 | + :type p: :class:`numpy.ndarray` |
| 871 | + :return: The shape function values at `p` [1 x 6]. |
| 872 | + :rtype: :class:`numpy.ndarray` |
| 873 | + """ |
| 874 | + eta = p[0] |
| 875 | + xi = p[1] |
| 876 | + zeta = p[2] |
| 877 | + |
| 878 | + return np.array( |
| 879 | + [ |
| 880 | + eta * (2 * eta - 1), |
| 881 | + xi * (2 * xi - 1), |
| 882 | + zeta * (2 * zeta - 1), |
| 883 | + 4 * eta * xi, |
| 884 | + 4 * xi * zeta, |
| 885 | + 4 * eta * zeta, |
| 886 | + ] |
| 887 | + ) |
| 888 | + |
| 889 | + |
678 | 890 | def extrapolate_to_nodes(w): |
679 | 891 | """Extrapolates results at six Gauss points to the six nodes of a quadratic triangular element. |
680 | 892 |
|
|
0 commit comments