PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
gmsh_particles Namespace Reference

Functions

 circle_mesh_symmetric (xc=[0., 0., 0.], r=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 ellipse_mesh_symmetric (xc=[0., 0., 0.], rx=1., ry=0.5, h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 sphere_mesh_symmetric (xc=[0., 0., 0.], r=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 rectangle_mesh_symmetric (xc=[0., 0., 0.], Lx=1., Ly=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 hexagon_mesh_symmetric (xc=[0., 0., 0.], r=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 drum2d_mesh_symmetric (xc=[0., 0., 0.], r=1., width=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 polygon_mesh_symmetric (points, theta, xc=[0., 0., 0.], h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 cylindrical2d_wall_mesh (center=[0., 0., 0.], outer_radius=1.0, inner_radius=0.8, bar_width=0.2, bar_length=0.3, h=0.1, filename='mesh', vtk_out=False)
 
 triangle_mesh_symmetric (xc=[0., 0., 0.], r=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 cuboid_mesh_symmetric (xc=[0., 0., 0.], Lx=1., Ly=1., Lz=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 ellipsoid_mesh_symmetric (xc=[0., 0., 0.], rx=1., ry=0.5, rz=0.3, h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 cylinder_mesh_symmetric (xc=[0., 0., 0.], r=1., h=1., mesh_size=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 annulus_circle_mesh_symmetric (xc=[0., 0., 0.], r_outer=1., r_inner=0.5, h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 annulus_rectangle_mesh (xc=[0., 0., 0.], Lx=1., Ly=1., hole_Lx=0.3, hole_Ly=0.3, h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True)
 
 open_rectangle_mesh (xc=[0., 0., 0.], Lx=1., Ly=1., hole_Lx=0.5, hole_Ly=0.5, h=0.1, filename='mesh', vtk_out=False)
 
 open_pipe_mesh (xc=[0., 0., 0.], axis=[0., 0., 1.], length=2., outer_radius=1., wall_thickness=0.1, h=0.1, filename='mesh', vtk_out=False)
 

Variables

str inp_dir = './'
 
list test_meshes = ['circle', 'ellipse', 'sphere', 'cuboid', 'ellipsoid', 'rectangle', 'hexagon', 'drum2d', 'triangle', 'polygon', 'cylindrical2d_wall', 'cylinder', 'annulus_circle', 'annulus_rectangle', 'open_rectangle', 'open_pipe']
 
int symm_flag = 0
 
 xc
 
 r
 
 h
 
 filename
 
 vtk_out
 
 True
 
 symmetric_mesh
 
 rx
 
 ry
 
 Lx
 
 Ly
 
 width
 
int theta = np.pi/6.
 
 R
 
 a
 
list v1 = [0., 0., 0.]
 
list v2 = [R*np.cos(0.5*theta), -R*np.sin(0.5*theta), 0.]
 
list v4 = [R*np.cos(0.5*theta), R*np.sin(0.5*theta), 0.]
 
list v3 = [R + a, 0., 0.]
 
 points
 
 center
 
 outer_radius
 
 inner_radius
 
 bar_width
 
 bar_length
 
 Lz
 
 rz
 
 mesh_size
 
 r_outer
 
 r_inner
 
 hole_Lx
 
 hole_Ly
 
 mesh
 
 axis
 
 length
 
 wall_thickness
 

Function Documentation

◆ annulus_circle_mesh_symmetric()

gmsh_particles.annulus_circle_mesh_symmetric (   xc = [0., 0., 0.],
  r_outer = 1.,
  r_inner = 0.5,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create an annulus (ring) mesh, either symmetric or full.
xc - center point coordinates [x, y, z]
r_outer - radius of outer circle
r_inner - radius of inner circle
h - mesh size
symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full annulus

Definition at line 989 of file gmsh_particles.py.

989def annulus_circle_mesh_symmetric(xc = [0., 0., 0.], r_outer = 1., r_inner = 0.5, h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
990 """
991 Create an annulus (ring) mesh, either symmetric or full.
992 xc - center point coordinates [x, y, z]
993 r_outer - radius of outer circle
994 r_inner - radius of inner circle
995 h - mesh size
996 symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full annulus
997 """
998 if r_inner >= r_outer:
999 raise ValueError("Inner radius must be smaller than outer radius")
1000
1001 gmsh.initialize()
1002 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
1003
1004 if symmetric_mesh:
1005 xc_mesh = [0., 0., 0.]
1006
1007 # add points for the quadrant of annulus circle
1008 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h)
1009 p2 = gmsh.model.geo.addPoint(xc_mesh[0] + r_outer, xc_mesh[1], xc_mesh[2], h)
1010 p3 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + r_outer, xc_mesh[2], h)
1011 p4 = gmsh.model.geo.addPoint(xc_mesh[0] + r_inner, xc_mesh[1], xc_mesh[2], h)
1012 p5 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + r_inner, xc_mesh[2], h)
1013
1014 l1 = gmsh.model.geo.addCircleArc(p2, p1, p3)
1015 l2 = gmsh.model.geo.addCircleArc(p4, p1, p5)
1016 l3 = gmsh.model.geo.addLine(p4, p2)
1017 l4 = gmsh.model.geo.addLine(p5, p3)
1018
1019 c1 = gmsh.model.geo.addCurveLoop([l3, l1, -l4, -l2])
1020 s1 = gmsh.model.geo.addPlaneSurface([c1])
1021
1022 gmsh.model.geo.synchronize()
1023 gmsh.model.mesh.generate(3)
1024
1025 m = get_gmsh_entities()
1026
1027 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # x-axis
1028 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # y-axis
1029 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # z-axis
1030
1031 # remove the duplicate nodes that will have been created on the internal
1032 # boundaries
1033 gmsh.model.mesh.removeDuplicateNodes()
1034
1035 # translate the mesh to the specified center coordinates
1036 gmsh_translate(xc)
1037
1038 else:
1039
1040 # Outer circle points
1041 p1 = gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h)
1042 p2 = gmsh.model.geo.addPoint(xc[0] + r_outer, xc[1], xc[2], h)
1043 p3 = gmsh.model.geo.addPoint(xc[0], xc[1] + r_outer, xc[2], h)
1044 p4 = gmsh.model.geo.addPoint(xc[0] - r_outer, xc[1], xc[2], h)
1045 p5 = gmsh.model.geo.addPoint(xc[0], xc[1] - r_outer, xc[2], h)
1046
1047 # Inner circle points
1048 p6 = gmsh.model.geo.addPoint(xc[0] + r_inner, xc[1], xc[2], h)
1049 p7 = gmsh.model.geo.addPoint(xc[0], xc[1] + r_inner, xc[2], h)
1050 p8 = gmsh.model.geo.addPoint(xc[0] - r_inner, xc[1], xc[2], h)
1051 p9 = gmsh.model.geo.addPoint(xc[0], xc[1] - r_inner, xc[2], h)
1052
1053 # Create circular arcs for outer circle
1054 c1 = gmsh.model.geo.addCircleArc(p2, p1, p3)
1055 c2 = gmsh.model.geo.addCircleArc(p3, p1, p4)
1056 c3 = gmsh.model.geo.addCircleArc(p4, p1, p5)
1057 c4 = gmsh.model.geo.addCircleArc(p5, p1, p2)
1058
1059 # Create circular arcs for inner circle
1060 c5 = gmsh.model.geo.addCircleArc(p6, p1, p7)
1061 c6 = gmsh.model.geo.addCircleArc(p7, p1, p8)
1062 c7 = gmsh.model.geo.addCircleArc(p8, p1, p9)
1063 c8 = gmsh.model.geo.addCircleArc(p9, p1, p6)
1064
1065 # Create curve loops
1066 outer_loop = gmsh.model.geo.addCurveLoop([c1, c2, c3, c4])
1067 inner_loop = gmsh.model.geo.addCurveLoop([c5, c6, c7, c8])
1068
1069 #gmsh.fltk.run()
1070
1071 # surface
1072 s1 = gmsh.model.geo.addPlaneSurface([outer_loop, inner_loop])
1073
1074 # designate surface as physical surface
1075 gmsh.model.addPhysicalGroup(2, [s1], 1)
1076
1077 # Synchronize and generate mesh
1078 gmsh.model.geo.synchronize()
1079 gmsh.model.mesh.generate(2)
1080
1081 # Write output files
1082 gmsh.write(filename + '.msh')
1083 if vtk_out:
1084 gmsh.write(filename + '.vtk')
1085
1086 gmsh.finalize()
1087

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ annulus_rectangle_mesh()

gmsh_particles.annulus_rectangle_mesh (   xc = [0., 0., 0.],
  Lx = 1.,
  Ly = 1.,
  hole_Lx = 0.3,
  hole_Ly = 0.3,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a rectangular mesh with a rectangular hole in the center (annulus rectangle).

Parameters
----------
xc : list
    Center coordinates [x, y, z]
Lx : float
    Length of outer rectangle in x direction
Ly : float
    Length of outer rectangle in y direction
hole_Lx : float
    Length of inner hole in x direction
hole_Ly : float
    Length of inner hole in y direction
h : float
    Mesh size
filename : str
    Output filename without extension
vtk_out : bool
    If True, also writes a VTK file
symmetric_mesh : bool
    If True, creates 1/4 mesh and mirrors it. If False, creates full mesh

Definition at line 1088 of file gmsh_particles.py.

1088def annulus_rectangle_mesh(xc=[0., 0., 0.], Lx=1., Ly=1., hole_Lx=0.3, hole_Ly=0.3, h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True):
1089 """
1090 Create a rectangular mesh with a rectangular hole in the center (annulus rectangle).
1091
1092 Parameters
1093 ----------
1094 xc : list
1095 Center coordinates [x, y, z]
1096 Lx : float
1097 Length of outer rectangle in x direction
1098 Ly : float
1099 Length of outer rectangle in y direction
1100 hole_Lx : float
1101 Length of inner hole in x direction
1102 hole_Ly : float
1103 Length of inner hole in y direction
1104 h : float
1105 Mesh size
1106 filename : str
1107 Output filename without extension
1108 vtk_out : bool
1109 If True, also writes a VTK file
1110 symmetric_mesh : bool
1111 If True, creates 1/4 mesh and mirrors it. If False, creates full mesh
1112 """
1113 # Initialize gmsh
1114 gmsh.initialize()
1115 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
1116
1117 if symmetric_mesh:
1118 # Create 1/4 mesh at origin first, then mirror and translate
1119 xc_mesh = [0., 0., 0.]
1120
1121 # Create points for outer rectangle (1/4)
1122 p1 = gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*Lx, xc_mesh[1], xc_mesh[2], h)
1123 p2 = gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*Lx, xc_mesh[1] + 0.5*Ly, xc_mesh[2], h)
1124 p3 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + 0.5*Ly, xc_mesh[2], h)
1125 p4 = gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*hole_Lx, xc_mesh[1], xc_mesh[2], h)
1126 p5 = gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*hole_Lx, xc_mesh[1] + 0.5*hole_Ly, xc_mesh[2], h)
1127 p6 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + 0.5*hole_Ly, xc_mesh[2], h)
1128
1129 # Create lines
1130 l1 = gmsh.model.geo.addLine(p1, p2)
1131 l2 = gmsh.model.geo.addLine(p2, p3)
1132 l3 = gmsh.model.geo.addLine(p3, p6)
1133 l4 = gmsh.model.geo.addLine(p6, p5)
1134 l5 = gmsh.model.geo.addLine(p5, p4)
1135 l6 = gmsh.model.geo.addLine(p4, p1)
1136
1137 # Create curve loops
1138 cl = gmsh.model.geo.addCurveLoop([l1, l2, l3, l4, l5, l6])
1139
1140 # Create plane surface with hole
1141 s = gmsh.model.geo.addPlaneSurface([cl])
1142
1143 # Synchronize and generate mesh
1144 gmsh.model.geo.synchronize()
1145 gmsh.model.mesh.generate(2)
1146
1147 # Get mesh data for mirroring
1148 m = get_gmsh_entities()
1149
1150 # Mirror the mesh in all quadrants
1151 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # Mirror across y-axis
1152 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # Mirror across x-axis
1153 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # Mirror across origin
1154
1155 # Remove duplicate nodes
1156 gmsh.model.mesh.removeDuplicateNodes()
1157
1158 # Translate to specified center coordinates
1159 gmsh_translate(xc)
1160
1161 else:
1162 # Create points for outer rectangle
1163 p1 = gmsh.model.geo.addPoint(xc[0] - 0.5*Lx, xc[1] - 0.5*Ly, xc[2], h) # Bottom left
1164 p2 = gmsh.model.geo.addPoint(xc[0] + 0.5*Lx, xc[1] - 0.5*Ly, xc[2], h) # Bottom right
1165 p3 = gmsh.model.geo.addPoint(xc[0] + 0.5*Lx, xc[1] + 0.5*Ly, xc[2], h) # Top right
1166 p4 = gmsh.model.geo.addPoint(xc[0] - 0.5*Lx, xc[1] + 0.5*Ly, xc[2], h) # Top left
1167
1168 # Create points for inner hole
1169 p5 = gmsh.model.geo.addPoint(xc[0] - 0.5*hole_Lx, xc[1] - 0.5*hole_Ly, xc[2], h) # Bottom left
1170 p6 = gmsh.model.geo.addPoint(xc[0] + 0.5*hole_Lx, xc[1] - 0.5*hole_Ly, xc[2], h) # Bottom right
1171 p7 = gmsh.model.geo.addPoint(xc[0] + 0.5*hole_Lx, xc[1] + 0.5*hole_Ly, xc[2], h) # Top right
1172 p8 = gmsh.model.geo.addPoint(xc[0] - 0.5*hole_Lx, xc[1] + 0.5*hole_Ly, xc[2], h) # Top left
1173
1174 # Create lines for outer rectangle
1175 l1 = gmsh.model.geo.addLine(p1, p2) # Bottom
1176 l2 = gmsh.model.geo.addLine(p2, p3) # Right
1177 l3 = gmsh.model.geo.addLine(p3, p4) # Top
1178 l4 = gmsh.model.geo.addLine(p4, p1) # Left
1179
1180 # Create lines for inner hole
1181 l5 = gmsh.model.geo.addLine(p5, p6) # Bottom
1182 l6 = gmsh.model.geo.addLine(p6, p7) # Right
1183 l7 = gmsh.model.geo.addLine(p7, p8) # Top
1184 l8 = gmsh.model.geo.addLine(p8, p5) # Left
1185
1186 # Create curve loops
1187 outer_loop = gmsh.model.geo.addCurveLoop([l1, l2, l3, l4])
1188 inner_loop = gmsh.model.geo.addCurveLoop([l5, l6, l7, l8])
1189
1190 # Create plane surface with hole
1191 s1 = gmsh.model.geo.addPlaneSurface([outer_loop, inner_loop])
1192
1193 # Synchronize and generate mesh
1194 gmsh.model.geo.synchronize()
1195 gmsh.model.mesh.generate(2)
1196
1197 # Write mesh files
1198 gmsh.write(filename + '.msh')
1199 if vtk_out:
1200 gmsh.write(filename + '.vtk')
1201
1202 # Finalize gmsh
1203 gmsh.finalize()
1204

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ circle_mesh_symmetric()

gmsh_particles.circle_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a circular mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
r - radius of the circle
h - mesh size
symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full circle

Definition at line 11 of file gmsh_particles.py.

11def circle_mesh_symmetric(xc = [0., 0., 0.], r = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
12 """
13 Create a circular mesh, either symmetric or full.
14 xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
15 r - radius of the circle
16 h - mesh size
17 symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full circle
18 """
19
20 gmsh.initialize()
21 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
22
23 if symmetric_mesh:
24 # we first assume circle center is at origin and then translate the symmetric mesh to the desired location
25 xc_mesh = [0., 0., 0.]
26 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h)
27 p2 = gmsh.model.geo.addPoint(xc_mesh[0] + r, xc_mesh[1], xc_mesh[2], h)
28 p3 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + r, xc_mesh[2], h)
29
30 l1 = gmsh.model.geo.addCircleArc(p2, p1, p3)
31 l2 = gmsh.model.geo.addLine(p1, p2)
32 l3 = gmsh.model.geo.addLine(p3, p1)
33
34 c1 = gmsh.model.geo.addCurveLoop([l2, l1, l3])
35
36 p1 = gmsh.model.geo.addPlaneSurface([c1])
37
38 #gmsh.model.occ.addBox(0,0,0, 1,0.5,0.5)
39 gmsh.model.geo.synchronize()
40 # gmsh.model.mesh.setSize(gmsh.model.getEntities(0), 0.1)
41 # gmsh.model.mesh.setSize([(0, 2)], 0.01)
42 gmsh.model.mesh.generate(3)
43 # gmsh.write('mesh.vtk')
44
45 # get the mesh data
46 m = get_gmsh_entities()
47
48 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # x-axis
49 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # y-axis
50 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # z-axis
51
52 # remove the duplicate nodes that will have been created on the internal
53 # boundaries
54 gmsh.model.mesh.removeDuplicateNodes()
55
56 # translate the mesh to the specified center coordinates
57 gmsh_translate(xc)
58
59 else:
60 # Create circle using built-in gmsh circle
61 c = gmsh.model.occ.addCircle(xc[0], xc[1], xc[2], r)
62
63 # Create curve loop from circle
64 cl = gmsh.model.occ.addCurveLoop([c])
65
66 # Create surface from curve loop
67 s = gmsh.model.occ.addPlaneSurface([cl])
68
69 # Add center point
70 p = gmsh.model.occ.addPoint(xc[0], xc[1], xc[2], h)
71
72 # Synchronize
73 gmsh.model.occ.synchronize()
74
75 # Embed center point in surface
76 gmsh.model.mesh.embed(0, [p], 2, s)
77
78 # generate mesh
79 gmsh.model.mesh.generate(3)
80
81 # write
82 gmsh.write(filename + '.msh')
83 if vtk_out:
84 gmsh.write(filename + '.vtk')
85
86 gmsh.finalize()
87

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Referenced by circle_wall_input_using_symmetric_mesh.create_input_file().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cuboid_mesh_symmetric()

gmsh_particles.cuboid_mesh_symmetric (   xc = [0., 0., 0.],
  Lx = 1.,
  Ly = 1.,
  Lz = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a cuboid mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left-back corner for non-symmetric mesh
Lx - length in x direction
Ly - length in y direction
Lz - length in z direction
h - mesh size
symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full cuboid

Definition at line 739 of file gmsh_particles.py.

739def cuboid_mesh_symmetric(xc = [0., 0., 0.], Lx = 1., Ly = 1., Lz = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
740 """
741 Create a cuboid mesh, either symmetric or full.
742 xc - center point for symmetric mesh, or bottom-left-back corner for non-symmetric mesh
743 Lx - length in x direction
744 Ly - length in y direction
745 Lz - length in z direction
746 h - mesh size
747 symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full cuboid
748 """
749 gmsh.initialize()
750 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
751
752 if symmetric_mesh:
753 # Create 1/8th of cuboid at origin, then mirror it
754 xc_mesh = [0., 0., 0.]
755
756 # Create box for first octant
757 box = gmsh.model.occ.addBox(xc_mesh[0], xc_mesh[1], xc_mesh[2],
758 0.5*Lx, 0.5*Ly, 0.5*Lz)
759
760 gmsh.model.occ.synchronize()
761
762 # Set mesh size
763 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
764
765 # Generate 3D mesh
766 gmsh.model.mesh.generate(3)
767
768 # Get mesh data for mirroring
769 m = get_gmsh_entities()
770
771 # Mirror the mesh to create full cuboid
772 # First create the other octants in the positive z hemisphere
773 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # Mirror across yz plane
774 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # Mirror across xz plane
775 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # Mirror across z axis
776
777 # Then mirror the entire upper hemisphere to create lower hemisphere
778 gmsh_transform(m, 4000, 4000000, 4000000, 1, 1, -1) # Mirror across xy plane
779 gmsh_transform(m, 5000, 5000000, 5000000, -1, 1, -1) # Mirror across yz plane
780 gmsh_transform(m, 6000, 6000000, 6000000, 1, -1, -1) # Mirror across xz plane
781 gmsh_transform(m, 7000, 7000000, 7000000, -1, -1, -1) # Mirror across z axis
782
783 # Remove duplicate nodes
784 gmsh.model.mesh.removeDuplicateNodes()
785
786 # Translate to specified center coordinates
787 gmsh_translate(xc)
788
789 else:
790 # Create full cuboid directly using OpenCASCADE
791 # Create box centered at xc
792 box = gmsh.model.occ.addBox(xc[0] - 0.5*Lx, xc[1] - 0.5*Ly, xc[2] - 0.5*Lz,
793 Lx, Ly, Lz)
794
795 # Add center point
796 p1 = gmsh.model.occ.addPoint(xc[0], xc[1], xc[2], h)
797
798 # Synchronize
799 gmsh.model.occ.synchronize()
800
801 # Embed center point in volume
802 gmsh.model.mesh.embed(0, [p1], 3, box)
803
804 # Set mesh size
805 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
806
807 # Generate 3D mesh
808 gmsh.model.mesh.generate(3)
809
810 # Write output files
811 gmsh.write(filename + '.msh')
812 if vtk_out:
813 gmsh.write(filename + '.vtk')
814
815 gmsh.finalize()
816

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ cylinder_mesh_symmetric()

gmsh_particles.cylinder_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  h = 1.,
  mesh_size = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a cylinder mesh, either symmetric or full.
xc - center point coordinates [x, y, z] (center of bottom face)
r - radius of cylinder
h - height of cylinder
mesh_size - mesh element size
symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full cylinder

Definition at line 901 of file gmsh_particles.py.

901def cylinder_mesh_symmetric(xc = [0., 0., 0.], r = 1., h = 1., mesh_size = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
902 """
903 Create a cylinder mesh, either symmetric or full.
904 xc - center point coordinates [x, y, z] (center of bottom face)
905 r - radius of cylinder
906 h - height of cylinder
907 mesh_size - mesh element size
908 symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full cylinder
909 """
910 gmsh.initialize()
911 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
912
913 if symmetric_mesh:
914 # Create 1/8th of cylinder at origin, then mirror it
915 xc_mesh = [0., 0., 0.]
916
917 # Create full cylinder first
918 cyl = gmsh.model.occ.addCylinder(xc_mesh[0], xc_mesh[1], xc_mesh[2],
919 0, 0, h, # height vector in z direction
920 r)
921
922 # Create cutting boxes for first octant
923 cut_left = gmsh.model.occ.addBox(xc_mesh[0]-r, xc_mesh[1]-r, xc_mesh[2], r, 2*r, h)
924 cut_back = gmsh.model.occ.addBox(xc_mesh[0]-r, xc_mesh[1]-r, xc_mesh[2], 2*r, r, h)
925
926 # Cut cylinder to get first quarter
927 gmsh.model.occ.cut([(3,cyl)], [(3,cut_left), (3,cut_back)])
928
929 gmsh.model.occ.synchronize()
930
931 # Set mesh size
932 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), mesh_size)
933
934 # Generate 3D mesh
935 gmsh.model.mesh.generate(3)
936
937 # Get mesh data for mirroring
938 m = get_gmsh_entities()
939
940 # Mirror the mesh to create full cylinder
941 # First create the other quarters
942 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # Mirror across yz plane
943 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # Mirror across xz plane
944 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # Mirror across z axis
945
946 # Remove duplicate nodes
947 gmsh.model.mesh.removeDuplicateNodes()
948
949 # Translate to specified center coordinates
950 gmsh_translate(xc)
951
952 else:
953 # Create full cylinder directly using OpenCASCADE
954 cyl = gmsh.model.occ.addCylinder(xc[0], xc[1], xc[2], # base center
955 0, 0, h, # height vector in z direction
956 r) # radius
957
958 # Add center points at top and bottom faces
959 p1 = gmsh.model.occ.addPoint(xc[0], xc[1], xc[2], mesh_size) # bottom center
960 p2 = gmsh.model.occ.addPoint(xc[0], xc[1], xc[2] + h, mesh_size) # top center
961
962 # Synchronize
963 gmsh.model.occ.synchronize()
964
965 # Embed center points in respective faces
966 # First get all surfaces
967 surfaces = gmsh.model.getEntities(2)
968 # Find top and bottom surfaces (they are parallel to xy-plane)
969 for s in surfaces:
970 com = gmsh.model.occ.getCenterOfMass(s[0], s[1])
971 if abs(com[2] - xc[2]) < 1e-6: # bottom face
972 gmsh.model.mesh.embed(0, [p1], 2, s[1])
973 elif abs(com[2] - (xc[2] + h)) < 1e-6: # top face
974 gmsh.model.mesh.embed(0, [p2], 2, s[1])
975
976 # Set mesh size
977 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), mesh_size)
978
979 # Generate 3D mesh
980 gmsh.model.mesh.generate(3)
981
982 # Write output files
983 gmsh.write(filename + '.msh')
984 if vtk_out:
985 gmsh.write(filename + '.vtk')
986
987 gmsh.finalize()
988

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ cylindrical2d_wall_mesh()

gmsh_particles.cylindrical2d_wall_mesh (   center = [0., 0., 0.],
  outer_radius = 1.0,
  inner_radius = 0.8,
  bar_width = 0.2,
  bar_length = 0.3,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False 
)
Create a cylindrical wall mesh with bars.

Parameters:
center - Center point coordinates [x, y, z]
outer_radius - Outer radius of the wall
inner_radius - Inner radius of the wall
bar_width - Width of the bars
bar_length - Length of the bars
h - Mesh size
filename - Output filename
vtk_out - Whether to output VTK file

Definition at line 576 of file gmsh_particles.py.

577 bar_width=0.2, bar_length=0.3, h=0.1, filename='mesh', vtk_out=False):
578 """
579 Create a cylindrical wall mesh with bars.
580
581 Parameters:
582 center - Center point coordinates [x, y, z]
583 outer_radius - Outer radius of the wall
584 inner_radius - Inner radius of the wall
585 bar_width - Width of the bars
586 bar_length - Length of the bars
587 h - Mesh size
588 filename - Output filename
589 vtk_out - Whether to output VTK file
590 """
591 gmsh.initialize()
592 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
593
594 # Create points
595 # Center point
596 p1 = gmsh.model.geo.addPoint(center[0], center[1], center[2], h)
597
598 # Outer circle points
599 p2 = gmsh.model.geo.addPoint(center[0] + outer_radius, center[1], center[2], h)
600 p3 = gmsh.model.geo.addPoint(center[0], center[1] + outer_radius, center[2], h)
601 p4 = gmsh.model.geo.addPoint(center[0] - outer_radius, center[1], center[2], h)
602 p5 = gmsh.model.geo.addPoint(center[0], center[1] - outer_radius, center[2], h)
603
604 # Inner points with bar
605 p6 = gmsh.model.geo.addPoint(center[0] + inner_radius, center[1] + 0.5*bar_width, center[2], h)
606 p7 = gmsh.model.geo.addPoint(center[0] + inner_radius, center[1] - 0.5*bar_width, center[2], h)
607 p8 = gmsh.model.geo.addPoint(center[0], center[1] + inner_radius, center[2], h)
608 p9 = gmsh.model.geo.addPoint(center[0] - inner_radius, center[1], center[2], h)
609 p10 = gmsh.model.geo.addPoint(center[0], center[1] - inner_radius, center[2], h)
610
611 # Bar end points
612 p11 = gmsh.model.geo.addPoint(center[0] + inner_radius - bar_length, center[1] + 0.5*bar_width, center[2], h)
613 p12 = gmsh.model.geo.addPoint(center[0] + inner_radius - bar_length, center[1] - 0.5*bar_width, center[2], h)
614
615 # Create circular arcs for outer circle
616 c1 = gmsh.model.geo.addCircleArc(p2, p1, p3)
617 c2 = gmsh.model.geo.addCircleArc(p3, p1, p4)
618 c3 = gmsh.model.geo.addCircleArc(p4, p1, p5)
619 c4 = gmsh.model.geo.addCircleArc(p5, p1, p2)
620
621 # Create circular arcs for inner circle
622 c5 = gmsh.model.geo.addCircleArc(p6, p1, p8)
623 c6 = gmsh.model.geo.addCircleArc(p8, p1, p9)
624 c7 = gmsh.model.geo.addCircleArc(p9, p1, p10)
625 c8 = gmsh.model.geo.addCircleArc(p10, p1, p7)
626
627 # Create lines for the bar
628 l1 = gmsh.model.geo.addLine(p6, p11)
629 l2 = gmsh.model.geo.addLine(p11, p12)
630 l3 = gmsh.model.geo.addLine(p12, p7)
631
632 # Create curve loops
633 outer_loop = gmsh.model.geo.addCurveLoop([c1, c2, c3, c4])
634 inner_loop = gmsh.model.geo.addCurveLoop([c5, c6, c7, c8, -l3, -l2, -l1])
635
636 # Create surface with hole
637 s1 = gmsh.model.geo.addPlaneSurface([outer_loop, inner_loop])
638
639 # Synchronize and generate mesh
640 gmsh.model.geo.synchronize()
641 gmsh.model.mesh.generate(2)
642
643 # Write output files
644 gmsh.write(filename + '.msh')
645 if vtk_out:
646 gmsh.write(filename + '.vtk')
647
648 gmsh.finalize()
649

◆ drum2d_mesh_symmetric()

gmsh_particles.drum2d_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  width = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a drum2d mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
r - radius of the drum2d
width - width of the drum2d
h - mesh size
symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full drum2d

Definition at line 419 of file gmsh_particles.py.

419def drum2d_mesh_symmetric(xc = [0., 0., 0.], r = 1., width = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
420 """
421 Create a drum2d mesh, either symmetric or full.
422 xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
423 r - radius of the drum2d
424 width - width of the drum2d
425 h - mesh size
426 symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full drum2d
427 """
428 gmsh.initialize()
429 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
430
431 if symmetric_mesh:
432 # Create 1/4 drum2d mesh at origin first, then rotate and translate
433 xc_mesh = [0., 0., 0.]
434 x_drum = get_ref_drum_points(xc_mesh, r, width)
435 v1, v2, v3 = x_drum[0], x_drum[1], x_drum[2]
436 v23 = [0.5*(v2[0] + v3[0]), 0.5*(v2[1] + v3[1]), 0.5*(v2[2] + v3[2])]
437
438 # 1/4th mesh is a polygon with 4 points (o, v1, v2, (v2+v3)/2)
439 # Create points for 1/4th mesh
440 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h) # Center point
441 p2 = gmsh.model.geo.addPoint(v1[0], v1[1], v1[2], h)
442 p3 = gmsh.model.geo.addPoint(v2[0], v2[1], v2[2], h)
443 p4 = gmsh.model.geo.addPoint(v23[0], v23[1], v23[2], h)
444
445 # Create lines
446 lines = []
447 points = [p1, p2, p3, p4, p1] # Add p1 again at end to complete loop
448 for i in range(4):
449 lines.append(gmsh.model.geo.addLine(points[i], points[i+1]))
450
451 # Create curve loop and surface
452 c1 = gmsh.model.geo.addCurveLoop(lines)
453 s1 = gmsh.model.geo.addPlaneSurface([c1])
454
455 gmsh.model.geo.synchronize()
456 gmsh.model.mesh.generate(3)
457
458 # Get mesh data for mirroring
459 m = get_gmsh_entities()
460
461 # Mirror mesh to create full drum2d
462 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # x-axis
463 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # y-axis
464 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # z-axis
465
466 # Remove duplicate nodes on internal boundaries
467 gmsh.model.mesh.removeDuplicateNodes()
468
469 # Translate to specified center
470 gmsh_translate(xc)
471
472 else:
473 # For non-symmetric mesh, create full drum2d directly
474 x_drum = get_ref_drum_points(xc, r, width)
475 # add points
476 points = []
477 points.append(gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h))
478 for i in range(len(x_drum)):
479 points.append(gmsh.model.geo.addPoint(x_drum[i][0], x_drum[i][1], x_drum[i][2], h))
480
481 # Create lines connecting vertices
482 lines = []
483 points.append(points[1]) # to close the loop
484 for i in range(len(x_drum)):
485 lines.append(gmsh.model.geo.addLine(points[i+1], points[i+2]))
486
487 # Create curve loop and surface
488 c1 = gmsh.model.geo.addCurveLoop(lines)
489 s1 = gmsh.model.geo.addPlaneSurface([c1])
490
491 # synchronize
492 gmsh.model.geo.synchronize()
493
494 # embed p1 (dim 0) in surface 1 (dim 2)
495 gmsh.model.mesh.embed(0, [points[0]], 2, s1)
496
497 # generate mesh
498 gmsh.model.mesh.generate(3)
499
500 # Write output files
501 gmsh.write(filename + '.msh')
502 if vtk_out:
503 gmsh.write(filename + '.vtk')
504
505 gmsh.finalize()
506

References geom_util.get_gmsh_entities(), geom_util.get_ref_drum_points(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ ellipse_mesh_symmetric()

gmsh_particles.ellipse_mesh_symmetric (   xc = [0., 0., 0.],
  rx = 1.,
  ry = 0.5,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create an elliptical mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
rx - radius in x direction
ry - radius in y direction
h - mesh size
symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full ellipse

Definition at line 88 of file gmsh_particles.py.

88def ellipse_mesh_symmetric(xc = [0., 0., 0.], rx = 1., ry = 0.5, h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
89 """
90 Create an elliptical mesh, either symmetric or full.
91 xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
92 rx - radius in x direction
93 ry - radius in y direction
94 h - mesh size
95 symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full ellipse
96 """
97
98 gmsh.initialize()
99 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
100
101 if symmetric_mesh:
102 # we first assume ellipse center is at origin and then translate the symmetric mesh to the desired location
103 xc_mesh = [0., 0., 0.]
104 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h) # Center
105 p2 = gmsh.model.geo.addPoint(xc_mesh[0] + rx, xc_mesh[1], xc_mesh[2], h) # Right
106 p3 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] + ry, xc_mesh[2], h) # Top
107
108 # Create elliptical arc using 3 points
109 l1 = gmsh.model.geo.addEllipseArc(p2, p1, p2, p3) # Right to top
110 l2 = gmsh.model.geo.addLine(p1, p2) # Center to right
111 l3 = gmsh.model.geo.addLine(p3, p1) # Top to center
112
113 # Create curve loop and surface
114 c1 = gmsh.model.geo.addCurveLoop([l2, l1, l3])
115 s1 = gmsh.model.geo.addPlaneSurface([c1])
116
117 gmsh.model.geo.synchronize()
118 gmsh.model.mesh.generate(3)
119
120 # get the mesh data
121 m = get_gmsh_entities()
122
123 # Mirror the mesh to create full ellipse
124 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # x-axis
125 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # y-axis
126 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # z-axis
127
128 # Remove duplicate nodes on internal boundaries
129 gmsh.model.mesh.removeDuplicateNodes()
130
131 # Translate to specified center coordinates
132 gmsh_translate(xc)
133
134 else:
135 # Create points for full ellipse
136 p1 = gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h) # Center
137 p2 = gmsh.model.geo.addPoint(xc[0] + rx, xc[1], xc[2], h) # Right
138 p3 = gmsh.model.geo.addPoint(xc[0], xc[1] + ry, xc[2], h) # Top
139 p4 = gmsh.model.geo.addPoint(xc[0] - rx, xc[1], xc[2], h) # Left
140 p5 = gmsh.model.geo.addPoint(xc[0], xc[1] - ry, xc[2], h) # Bottom
141
142 # Create elliptical arcs
143 l1 = gmsh.model.geo.addEllipseArc(p2, p1, p2, p3) # Right to top
144 l2 = gmsh.model.geo.addEllipseArc(p3, p1, p3, p4) # Top to left
145 l3 = gmsh.model.geo.addEllipseArc(p4, p1, p4, p5) # Left to bottom
146 l4 = gmsh.model.geo.addEllipseArc(p5, p1, p5, p2) # Bottom to right
147
148 # Create curve loop and surface
149 c1 = gmsh.model.geo.addCurveLoop([l1, l2, l3, l4])
150 s1 = gmsh.model.geo.addPlaneSurface([c1])
151
152 # synchronize
153 gmsh.model.geo.synchronize()
154
155 # embed p0 (dim 0) in surface 1 (dim 2)
156 gmsh.model.mesh.embed(0, [p1], 2, s1)
157
158 # generate mesh
159 gmsh.model.mesh.generate(3)
160
161 # Write output files
162 gmsh.write(filename + '.msh')
163 if vtk_out:
164 gmsh.write(filename + '.vtk')
165
166 gmsh.finalize()
167

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ ellipsoid_mesh_symmetric()

gmsh_particles.ellipsoid_mesh_symmetric (   xc = [0., 0., 0.],
  rx = 1.,
  ry = 0.5,
  rz = 0.3,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create an ellipsoid mesh, either symmetric or full.
xc - center point coordinates [x, y, z]
rx - radius in x direction
ry - radius in y direction
rz - radius in z direction
h - mesh size
symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full ellipsoid

Definition at line 817 of file gmsh_particles.py.

817def ellipsoid_mesh_symmetric(xc = [0., 0., 0.], rx = 1., ry = 0.5, rz = 0.3, h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
818 """
819 Create an ellipsoid mesh, either symmetric or full.
820 xc - center point coordinates [x, y, z]
821 rx - radius in x direction
822 ry - radius in y direction
823 rz - radius in z direction
824 h - mesh size
825 symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full ellipsoid
826 """
827 gmsh.initialize()
828 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
829
830 if symmetric_mesh:
831 # Create 1/8th of ellipsoid at origin, then mirror it
832 xc_mesh = [0., 0., 0.]
833
834 # Create sphere of radius 1
835 sphere = gmsh.model.occ.addSphere(xc_mesh[0], xc_mesh[1], xc_mesh[2], 1.0)
836
837 # Scale it to create ellipsoid
838 gmsh.model.occ.dilate([(3, sphere)], xc_mesh[0], xc_mesh[1], xc_mesh[2], rx, ry, rz)
839
840 # Create cutting boxes for first octant
841 cut_left = gmsh.model.occ.addBox(xc_mesh[0]-rx, xc_mesh[1]-ry, xc_mesh[2]-rz, rx, 2*ry, 2*rz)
842 cut_back = gmsh.model.occ.addBox(xc_mesh[0]-rx, xc_mesh[1]-ry, xc_mesh[2]-rz, 2*rx, ry, 2*rz)
843 cut_bottom = gmsh.model.occ.addBox(xc_mesh[0]-rx, xc_mesh[1]-ry, xc_mesh[2]-rz, 2*rx, 2*ry, rz)
844
845 # Cut sphere to get first octant
846 gmsh.model.occ.cut([(3,sphere)], [(3,cut_left), (3,cut_back), (3,cut_bottom)])
847
848 gmsh.model.occ.synchronize()
849
850 # Set mesh size
851 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
852
853 # Generate 3D mesh
854 gmsh.model.mesh.generate(3)
855
856 # Get mesh data for mirroring
857 m = get_gmsh_entities()
858
859 # Mirror the mesh to create full ellipsoid
860 # First create the other octants in the positive z hemisphere
861 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1) # Mirror across yz plane
862 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1) # Mirror across xz plane
863 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1) # Mirror across z axis
864
865 # Then mirror the entire upper hemisphere to create lower hemisphere
866 gmsh_transform(m, 4000, 4000000, 4000000, 1, 1, -1) # Mirror across xy plane
867 gmsh_transform(m, 5000, 5000000, 5000000, -1, 1, -1) # Mirror across yz plane
868 gmsh_transform(m, 6000, 6000000, 6000000, 1, -1, -1) # Mirror across xz plane
869 gmsh_transform(m, 7000, 7000000, 7000000, -1, -1, -1) # Mirror across z axis
870
871 # Remove duplicate nodes
872 gmsh.model.mesh.removeDuplicateNodes()
873
874 # Translate to specified center coordinates
875 gmsh_translate(xc)
876
877 else:
878 # Create full ellipsoid directly using OpenCASCADE
879 # First create sphere of radius 1
880 sphere = gmsh.model.occ.addSphere(xc[0], xc[1], xc[2], 1.0)
881
882 # Scale it to create ellipsoid
883 gmsh.model.occ.dilate([(3, sphere)], xc[0], xc[1], xc[2], rx, ry, rz)
884
885 # Synchronize
886 gmsh.model.occ.synchronize()
887
888 # Set mesh size
889 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
890
891 # Generate 3D mesh
892 gmsh.model.mesh.generate(3)
893
894 # Write output files
895 gmsh.write(filename + '.msh')
896 if vtk_out:
897 gmsh.write(filename + '.vtk')
898
899 gmsh.finalize()
900

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ hexagon_mesh_symmetric()

gmsh_particles.hexagon_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a hexagon mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
r - radius of the hexagon
h - mesh size
symmetric_mesh - if True, creates 1/6 mesh and mirrors it. If False, creates full hexagon

Definition at line 338 of file gmsh_particles.py.

338def hexagon_mesh_symmetric(xc = [0., 0., 0.], r = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
339 """
340 Create a hexagon mesh, either symmetric or full.
341 xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
342 r - radius of the hexagon
343 h - mesh size
344 symmetric_mesh - if True, creates 1/6 mesh and mirrors it. If False, creates full hexagon
345 """
346 gmsh.initialize()
347 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
348
349 if symmetric_mesh:
350 # Create 1/6 hexagon mesh at origin first, then rotate and translate
351 xc_mesh = [0., 0., 0.]
352 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h) # Center (origin)
353 p2 = gmsh.model.geo.addPoint(xc_mesh[0] + r, xc_mesh[1], xc_mesh[2], h) # Point on x-axis (p1)
354 p3 = gmsh.model.geo.addPoint(xc_mesh[0] + r*np.cos(np.pi/3), xc_mesh[1] + r*np.sin(np.pi/3), xc_mesh[2], h) # p2
355
356 # Create lines for 1/6th of hexagon
357 l1 = gmsh.model.geo.addLine(p1, p2) # Center to p1
358 l2 = gmsh.model.geo.addLine(p2, p3) # p1 to p2
359 l3 = gmsh.model.geo.addLine(p3, p1) # p2 to center
360
361 # Create curve loop and surface
362 c1 = gmsh.model.geo.addCurveLoop([l1, l2, l3])
363 s1 = gmsh.model.geo.addPlaneSurface([c1])
364
365 gmsh.model.geo.synchronize()
366 gmsh.model.mesh.generate(3)
367
368 # Get mesh data for rotation
369 m = get_gmsh_entities()
370
371 # Rotate the mesh 5 times by 60 degrees
372 for i in range(5):
373 angle = (i + 1) * np.pi/3
374 gmsh_transform_general(m, 1000*(i+1), 1000000*(i+1), 1000000*(i+1), angle)
375
376 # Remove duplicate nodes
377 gmsh.model.mesh.removeDuplicateNodes()
378
379 # Translate to specified center
380 gmsh_translate(xc)
381
382 else:
383 # For non-symmetric mesh, create full hexagon directly
384 points = []
385 for i in range(6):
386 angle = i * np.pi/3
387 px = xc[0] + r * np.cos(angle)
388 py = xc[1] + r * np.sin(angle)
389 points.append(gmsh.model.geo.addPoint(px, py, xc[2], h))
390
391 points.append(gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h))
392
393 # Create lines connecting points
394 lines = []
395 for i in range(5):
396 lines.append(gmsh.model.geo.addLine(points[i], points[i+1]))
397 lines.append(gmsh.model.geo.addLine(points[5], points[0]))
398
399 # Create curve loop and surface
400 c1 = gmsh.model.geo.addCurveLoop(lines)
401 s1 = gmsh.model.geo.addPlaneSurface([c1])
402
403 # synchronize
404 gmsh.model.geo.synchronize()
405
406 # embed p1 (dim 0) in surface 1 (dim 2)
407 gmsh.model.mesh.embed(0, [points[-1]], 2, s1)
408
409 # generate mesh
410 gmsh.model.mesh.generate(3)
411
412 # Write output files
413 gmsh.write(filename + '.msh')
414 if vtk_out:
415 gmsh.write(filename + '.vtk')
416
417 gmsh.finalize()
418

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform_general(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ open_pipe_mesh()

gmsh_particles.open_pipe_mesh (   xc = [0., 0., 0.],
  axis = [0., 0., 1.],
  length = 2.,
  outer_radius = 1.,
  wall_thickness = 0.1,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False 
)
Create a 3D pipe mesh with specified axis, closed bottom and open top.
The pipe is created by:
1. Creating an outer cylinder
2. Subtracting an inner cylinder to create walls
3. Subtracting the top surface to create the opening

Parameters
----------
xc : list
    Center coordinates [x, y, z] of the base center
axis : list
    Axis vector defining pipe orientation [ax, ay, az]
length : float
    Length of the pipe along axis
outer_radius : float
    Outer radius of the pipe
wall_thickness : float
    Thickness of the pipe wall and bottom
h : float
    Mesh size
filename : str
    Output filename without extension
vtk_out : bool
    If True, also writes a VTK file

Definition at line 1266 of file gmsh_particles.py.

1266def open_pipe_mesh(xc=[0., 0., 0.], axis=[0., 0., 1.], length=2., outer_radius=1., wall_thickness=0.1, h=0.1, filename='mesh', vtk_out=False):
1267 """
1268 Create a 3D pipe mesh with specified axis, closed bottom and open top.
1269 The pipe is created by:
1270 1. Creating an outer cylinder
1271 2. Subtracting an inner cylinder to create walls
1272 3. Subtracting the top surface to create the opening
1273
1274 Parameters
1275 ----------
1276 xc : list
1277 Center coordinates [x, y, z] of the base center
1278 axis : list
1279 Axis vector defining pipe orientation [ax, ay, az]
1280 length : float
1281 Length of the pipe along axis
1282 outer_radius : float
1283 Outer radius of the pipe
1284 wall_thickness : float
1285 Thickness of the pipe wall and bottom
1286 h : float
1287 Mesh size
1288 filename : str
1289 Output filename without extension
1290 vtk_out : bool
1291 If True, also writes a VTK file
1292 """
1293 gmsh.initialize()
1294 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
1295
1296 # Normalize axis vector
1297 axis_norm = np.sqrt(axis[0]**2 + axis[1]**2 + axis[2]**2)
1298 if axis_norm == 0:
1299 raise ValueError("Axis vector cannot be zero")
1300 axis = [x/axis_norm for x in axis]
1301
1302 # Create outer cylinder
1303 l = length + wall_thickness
1304 outer_cylinder = gmsh.model.occ.addCylinder(xc[0], xc[1], xc[2],
1305 axis[0]*l, axis[1]*l, axis[2]*l,
1306 outer_radius)
1307
1308 # Create inner cylinder (to be subtracted)
1309 inner_cylinder = gmsh.model.occ.addCylinder(xc[0]+wall_thickness*axis[0], xc[1]+wall_thickness*axis[1], xc[2]+wall_thickness*axis[2],
1310 axis[0]*l, axis[1]*l, axis[2]*l,
1311 outer_radius - wall_thickness)
1312
1313 # Create a box slightly larger than the cylinder at the top to cut off the top surface
1314 # First create a box aligned with coordinate axes
1315 dx = 2 * outer_radius
1316 box = gmsh.model.occ.addBox(xc[0] - dx, xc[1] - dx, xc[2] + l - wall_thickness,
1317 2*dx, 2*dx, 2*wall_thickness)
1318
1319 # If axis is not aligned with z, rotate the box to align with axis
1320 if not (abs(axis[0]) < 1e-10 and abs(axis[1]) < 1e-10):
1321 # Calculate rotation angle and axis
1322 z_axis = [0, 0, 1]
1323 rotation_axis = np.cross(z_axis, axis)
1324 rotation_angle = np.arccos(np.dot(z_axis, axis))
1325
1326 # Rotate the box around center point
1327 gmsh.model.occ.rotate([(3, box)],
1328 xc[0], xc[1], xc[2],
1329 rotation_axis[0], rotation_axis[1], rotation_axis[2],
1330 rotation_angle)
1331
1332 # First subtract inner cylinder from outer cylinder
1333 pipe = gmsh.model.occ.cut([(3, outer_cylinder)], [(3, inner_cylinder), (3, box)])
1334
1335 # Synchronize before meshing
1336 gmsh.model.occ.synchronize()
1337
1338 # Set mesh size
1339 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
1340
1341 # Generate 3D mesh
1342 gmsh.model.mesh.generate(3)
1343
1344 # Write mesh files
1345 gmsh.write(filename + '.msh')
1346 if vtk_out:
1347 gmsh.write(filename + '.vtk')
1348
1349 gmsh.finalize()
1350

◆ open_rectangle_mesh()

gmsh_particles.open_rectangle_mesh (   xc = [0., 0., 0.],
  Lx = 1.,
  Ly = 1.,
  hole_Lx = 0.5,
  hole_Ly = 0.5,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False 
)
Create a rectangular mesh with a rectangular hole in the center (annulus rectangle).

Parameters
----------
xc : list
    Center coordinates [x, y, z]
Lx : float
    Length of outer rectangle in x direction
Ly : float
    Length of outer rectangle in y direction
hole_Lx : float
    Length of inner hole in x direction
hole_Ly : float
    Length of inner hole in y direction
h : float
    Mesh size
filename : str
    Output filename without extension
vtk_out : bool
    If True, also writes a VTK file

Definition at line 1205 of file gmsh_particles.py.

1205def open_rectangle_mesh(xc = [0., 0., 0.], Lx = 1., Ly = 1., hole_Lx = 0.5, hole_Ly = 0.5, h = 0.1, filename = 'mesh', vtk_out = False):
1206 """
1207 Create a rectangular mesh with a rectangular hole in the center (annulus rectangle).
1208
1209 Parameters
1210 ----------
1211 xc : list
1212 Center coordinates [x, y, z]
1213 Lx : float
1214 Length of outer rectangle in x direction
1215 Ly : float
1216 Length of outer rectangle in y direction
1217 hole_Lx : float
1218 Length of inner hole in x direction
1219 hole_Ly : float
1220 Length of inner hole in y direction
1221 h : float
1222 Mesh size
1223 filename : str
1224 Output filename without extension
1225 vtk_out : bool
1226 If True, also writes a VTK file
1227 """
1228 # Initialize gmsh
1229 gmsh.initialize()
1230 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
1231
1232 # Annulus rectangle with top removed
1233 # xc is the center of the cocentric rectangles
1234 # we mesh the right half of the annulus rectangle
1235 xc_mesh = [0., 0., 0.]
1236 points = []
1237 points.append(gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*Lx, xc_mesh[1] - 0.5*Ly, xc_mesh[2], h))
1238 points.append(gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*Lx, xc_mesh[1] + 0.5*Ly, xc_mesh[2], h))
1239 points.append(gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*hole_Lx, xc_mesh[1] + 0.5*Ly, xc_mesh[2], h))
1240 points.append(gmsh.model.geo.addPoint(xc_mesh[0] + 0.5*hole_Lx, xc_mesh[1] - 0.5*hole_Ly, xc_mesh[2], h))
1241 points.append(gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] - 0.5*hole_Ly, xc_mesh[2], h))
1242 points.append(gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1] - 0.5*Ly, xc_mesh[2], h))
1243
1244 lines = []
1245 for i in range(len(points) - 1):
1246 lines.append(gmsh.model.geo.addLine(points[i], points[i + 1]))
1247 lines.append(gmsh.model.geo.addLine(points[-1], points[0]))
1248 cl = gmsh.model.geo.addCurveLoop(lines)
1249 s = gmsh.model.geo.addPlaneSurface([cl])
1250 gmsh.model.geo.synchronize()
1251 gmsh.model.mesh.generate(2)
1252
1253 m = get_gmsh_entities()
1254
1255 # mirror the mesh across the x-axis
1256 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1)
1257 gmsh.model.mesh.removeDuplicateNodes()
1258 gmsh_translate(xc)
1259
1260 gmsh.write(filename + '.msh')
1261 if vtk_out:
1262 gmsh.write(filename + '.vtk')
1263
1264 gmsh.finalize()
1265

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ polygon_mesh_symmetric()

gmsh_particles.polygon_mesh_symmetric (   points,
  theta,
  xc = [0., 0., 0.],
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a symmetric mesh by rotating a polygon segment n times, where n = 360/theta.

Parameters:
points - List of [x,y,z] coordinates defining the polygon segment. First point must be [0,0,0]
        and first two edges must form angle theta at origin
theta - Angle in radians between first two edges at origin. Must divide 2pi evenly.
xc - Center point for final translation
h - Mesh size
filename - Output filename
vtk_out - Whether to output VTK file

Definition at line 507 of file gmsh_particles.py.

507def polygon_mesh_symmetric(points, theta, xc=[0., 0., 0.], h=0.1, filename='mesh', vtk_out=False, symmetric_mesh = True):
508 """
509 Create a symmetric mesh by rotating a polygon segment n times, where n = 360/theta.
510
511 Parameters:
512 points - List of [x,y,z] coordinates defining the polygon segment. First point must be [0,0,0]
513 and first two edges must form angle theta at origin
514 theta - Angle in radians between first two edges at origin. Must divide 2pi evenly.
515 xc - Center point for final translation
516 h - Mesh size
517 filename - Output filename
518 vtk_out - Whether to output VTK file
519 """
520 # Validate inputs
521 if not points or len(points) < 3:
522 raise ValueError("Need at least 3 points to define a polygon segment")
523
524 if not np.allclose(points[0], [0., 0., 0.]):
525 raise ValueError("First point must be at origin [0,0,0]")
526
527 # Check if theta divides 360 evenly
528 n = int(round(2*np.pi/theta))
529 if not np.isclose(2*np.pi, n * theta):
530 raise ValueError(f"Angle {theta} degrees must divide 360 evenly")
531
532 gmsh.initialize()
533 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
534
535 # Create points for the polygon segment
536 gmsh_points = []
537 for pt in points:
538 gmsh_points.append(gmsh.model.geo.addPoint(pt[0], pt[1], pt[2], h))
539
540 # Create lines connecting points
541 lines = []
542 for i in range(len(gmsh_points)-1):
543 lines.append(gmsh.model.geo.addLine(gmsh_points[i], gmsh_points[i+1]))
544 # Close the polygon segment
545 lines.append(gmsh.model.geo.addLine(gmsh_points[-1], gmsh_points[0]))
546
547 # Create curve loop and surface
548 c1 = gmsh.model.geo.addCurveLoop(lines)
549 s1 = gmsh.model.geo.addPlaneSurface([c1])
550
551 # Synchronize and generate initial mesh
552 gmsh.model.geo.synchronize()
553 gmsh.model.mesh.generate(3)
554
555 # Get mesh data for rotation
556 m = get_gmsh_entities()
557
558 # Rotate n-1 times by theta
559 for i in range(1, n):
560 angle = i * theta # Convert to radians
561 gmsh_transform_general(m, 1000*i, 1000000*i, 1000000*i, angle)
562
563 # Remove duplicate nodes
564 gmsh.model.mesh.removeDuplicateNodes()
565
566 # Translate to final position if needed
567 gmsh_translate(xc)
568
569 # Write output files
570 gmsh.write(filename + '.msh')
571 if vtk_out:
572 gmsh.write(filename + '.vtk')
573
574 gmsh.finalize()
575

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform_general(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ rectangle_mesh_symmetric()

gmsh_particles.rectangle_mesh_symmetric (   xc = [0., 0., 0.],
  Lx = 1.,
  Ly = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a rectangular mesh, either symmetric or full.
xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
Lx - length in x direction 
Ly - length in y direction
h - mesh size
symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full rectangle

Definition at line 253 of file gmsh_particles.py.

253def rectangle_mesh_symmetric(xc = [0., 0., 0.], Lx = 1., Ly = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
254 """
255 Create a rectangular mesh, either symmetric or full.
256 xc - center point for symmetric mesh, or bottom-left corner for non-symmetric mesh
257 Lx - length in x direction
258 Ly - length in y direction
259 h - mesh size
260 symmetric_mesh - if True, creates 1/4 mesh and mirrors it. If False, creates full rectangle
261 """
262 # Initialize gmsh
263 gmsh.initialize()
264 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
265
266 if symmetric_mesh:
267 # center at origin
268 xc_mesh = [0., 0., 0.]
269
270 # Create points for 1/4 rectangle
271 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h)
272 p2 = gmsh.model.geo.addPoint(xc_mesh[0]+0.5*Lx, xc_mesh[1], xc_mesh[2], h)
273 p3 = gmsh.model.geo.addPoint(xc_mesh[0]+0.5*Lx, xc_mesh[1]+0.5*Ly, xc_mesh[2], h)
274 p4 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1]+0.5*Ly, xc_mesh[2], h)
275
276 # Create lines
277 l1 = gmsh.model.geo.addLine(p1, p2)
278 l2 = gmsh.model.geo.addLine(p2, p3)
279 l3 = gmsh.model.geo.addLine(p3, p4)
280 l4 = gmsh.model.geo.addLine(p4, p1)
281
282 # Create curve loop and surface
283 c1 = gmsh.model.geo.addCurveLoop([l1, l2, l3, l4])
284 s1 = gmsh.model.geo.addPlaneSurface([c1])
285
286 # Synchronize and generate mesh
287 gmsh.model.geo.synchronize()
288 gmsh.model.mesh.generate(3)
289
290 # Get the mesh data for mirroring
291 m = get_gmsh_entities()
292
293 # Mirror the mesh in all quadrants
294 gmsh_transform(m, 1000, 1000000, 1000000, -1, 1, 1)
295 gmsh_transform(m, 2000, 2000000, 2000000, 1, -1, 1)
296 gmsh_transform(m, 3000, 3000000, 3000000, -1, -1, 1)
297
298 # Remove duplicate nodes
299 gmsh.model.mesh.removeDuplicateNodes()
300
301 # translate the mesh to the specified center coordinates
302 gmsh_translate(xc)
303
304 else:
305 # Create points for full rectangle
306 p0 = gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h) # Center
307 p1 = gmsh.model.geo.addPoint(xc[0] - 0.5*Lx, xc[1] - 0.5*Ly, xc[2], h) # Bottom left
308 p2 = gmsh.model.geo.addPoint(xc[0] + 0.5*Lx, xc[1] - 0.5*Ly, xc[2], h) # Bottom right
309 p3 = gmsh.model.geo.addPoint(xc[0] + 0.5*Lx, xc[1] + 0.5*Ly, xc[2], h) # Top right
310 p4 = gmsh.model.geo.addPoint(xc[0] - 0.5*Lx, xc[1] + 0.5*Ly, xc[2], h) # Top left
311
312 # Create lines
313 l1 = gmsh.model.geo.addLine(p1, p2) # Bottom
314 l2 = gmsh.model.geo.addLine(p2, p3) # Right
315 l3 = gmsh.model.geo.addLine(p3, p4) # Top
316 l4 = gmsh.model.geo.addLine(p4, p1) # Left
317
318 # Create curve loop and surface
319 c1 = gmsh.model.geo.addCurveLoop([l1, l2, l3, l4])
320 s1 = gmsh.model.geo.addPlaneSurface([c1])
321
322 # synchronize
323 gmsh.model.geo.synchronize()
324
325 # embed p0 (dim 0) in surface 1 (dim 2)
326 gmsh.model.mesh.embed(0, [p0], 2, s1)
327
328 # generate mesh
329 gmsh.model.mesh.generate(3)
330
331 # write
332 gmsh.write(filename + '.msh')
333 if vtk_out:
334 gmsh.write(filename + '.vtk')
335
336 gmsh.finalize()
337

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Referenced by circle_wall_input_using_symmetric_mesh.create_input_file().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sphere_mesh_symmetric()

gmsh_particles.sphere_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a spherical mesh, either symmetric or full.
xc - center point coordinates [x, y, z]
r - radius of the sphere
h - mesh size
symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full sphere

Definition at line 168 of file gmsh_particles.py.

168def sphere_mesh_symmetric(xc = [0., 0., 0.], r = 1., h = 0.1, filename = 'mesh', vtk_out = False, symmetric_mesh = True):
169 """
170 Create a spherical mesh, either symmetric or full.
171 xc - center point coordinates [x, y, z]
172 r - radius of the sphere
173 h - mesh size
174 symmetric_mesh - if True, creates 1/8 mesh and mirrors it. If False, creates full sphere
175 """
176 gmsh.initialize()
177 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
178
179 if symmetric_mesh:
180 # Create 1/8th of sphere at origin, then mirror it
181 xc_mesh = [0., 0., 0.]
182
183 # Create first octant of sphere using OpenCASCADE
184 sphere = gmsh.model.occ.addSphere(xc_mesh[0], xc_mesh[1], xc_mesh[2], r)
185
186 # cubes to cut the sphere
187 cut_left = gmsh.model.occ.addBox(xc_mesh[0]-r, xc_mesh[1]-r, xc_mesh[2]-r, r, 2*r, 2*r)
188 cut_back = gmsh.model.occ.addBox(xc_mesh[0]-r, xc_mesh[1]-r, xc_mesh[2]-r, 2*r, r, 2*r)
189 cut_bottom = gmsh.model.occ.addBox(xc_mesh[0]-r, xc_mesh[1]-r, xc_mesh[2]-r, 2*r, 2*r, r)
190
191 # Cut sphere to get first octant using 3 planes
192 gmsh.model.occ.cut([(3,sphere)], [(3,cut_left), (3,cut_back), (3,cut_bottom)])
193
194 gmsh.model.occ.synchronize()
195
196 # Set mesh size
197 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
198
199 # Generate 3D mesh
200 gmsh.model.mesh.generate(3)
201
202 # gmsh.fltk.run()
203
204 # Get mesh data for mirroring
205 m = get_gmsh_entities()
206
207 # Mirror mesh in x direction
208 gmsh_transform(m, 1000, 10000, 100000, -1, 1, 1)
209 # Mirror in y direction
210 gmsh_transform(m, 2000, 20000, 200000, 1, -1, 1)
211 # Mirror in z direction
212 gmsh_transform(m, 3000, 30000, 300000, 1, 1, -1)
213 # Mirror in xy plane
214 gmsh_transform(m, 4000, 40000, 400000, -1, -1, 1)
215 # Mirror in yz plane
216 gmsh_transform(m, 5000, 50000, 500000, 1, -1, -1)
217 # Mirror in xz plane
218 gmsh_transform(m, 6000, 60000, 600000, -1, 1, -1)
219 # Mirror in xyz
220 gmsh_transform(m, 7000, 70000, 700000, -1, -1, -1)
221
222 # Translate to specified center if not at origin
223 gmsh_translate(xc)
224
225 else:
226 # Create points for full sphere
227 # Create full sphere using OpenCASCADE
228 sphere = gmsh.model.occ.addSphere(xc[0], xc[1], xc[2], r)
229
230 # Synchronize the OpenCASCADE geometry
231 gmsh.model.occ.synchronize()
232
233 # Set mesh size
234 gmsh.model.mesh.setSize(gmsh.model.getEntities(0), h)
235
236 # Synchronize
237 gmsh.model.geo.synchronize()
238
239 # Embed center point in volume
240 p1 = gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h)
241 gmsh.model.mesh.embed(0, [p1], 3, sphere)
242
243 # Generate 3D mesh
244 gmsh.model.mesh.generate(3)
245
246 # Write output files
247 gmsh.write(filename + '.msh')
248 if vtk_out:
249 gmsh.write(filename + '.vtk')
250
251 gmsh.finalize()
252

References geom_util.get_gmsh_entities(), geom_util.gmsh_transform(), and geom_util.gmsh_translate().

Here is the call graph for this function:

◆ triangle_mesh_symmetric()

gmsh_particles.triangle_mesh_symmetric (   xc = [0., 0., 0.],
  r = 1.,
  h = 0.1,
  filename = 'mesh',
  vtk_out = False,
  symmetric_mesh = True 
)
Create a triangle mesh, either symmetric or full.

Parameters:
xc - center point coordinates [x, y, z]
r - radius of the circumscribed circle
h - mesh size
filename - output filename
vtk_out - whether to output VTK file
symmetric_mesh - if True, creates 1/3 mesh and rotates it twice. If False, creates full triangle

Definition at line 650 of file gmsh_particles.py.

650def triangle_mesh_symmetric(xc=[0., 0., 0.], r=1., h=0.1, filename='mesh', vtk_out=False, symmetric_mesh=True):
651 """
652 Create a triangle mesh, either symmetric or full.
653
654 Parameters:
655 xc - center point coordinates [x, y, z]
656 r - radius of the circumscribed circle
657 h - mesh size
658 filename - output filename
659 vtk_out - whether to output VTK file
660 symmetric_mesh - if True, creates 1/3 mesh and rotates it twice. If False, creates full triangle
661 """
662 gmsh.initialize()
663 gmsh.option.setNumber("Mesh.MshFileVersion", 2.2)
664
665 if symmetric_mesh:
666 # Create 1/3 triangle mesh at origin first, then rotate and translate
667 xc_mesh = [0., 0., 0.]
668 # Get reference triangle points
669 x_tri = get_ref_triangle_points(xc_mesh, r)
670 v1, v2, v3 = x_tri[0], x_tri[1], x_tri[2]
671
672 # Create points for 1/3rd mesh (triangle from center to first two vertices)
673 p1 = gmsh.model.geo.addPoint(xc_mesh[0], xc_mesh[1], xc_mesh[2], h) # Center
674 p2 = gmsh.model.geo.addPoint(v1[0], v1[1], v1[2], h) # First vertex
675 p3 = gmsh.model.geo.addPoint(v2[0], v2[1], v2[2], h) # Second vertex
676
677 # Create lines
678 l1 = gmsh.model.geo.addLine(p1, p2) # Center to first vertex
679 l2 = gmsh.model.geo.addLine(p2, p3) # First to second vertex
680 l3 = gmsh.model.geo.addLine(p3, p1) # Second vertex back to center
681
682 # Create curve loop and surface
683 c1 = gmsh.model.geo.addCurveLoop([l1, l2, l3])
684 s1 = gmsh.model.geo.addPlaneSurface([c1])
685
686 gmsh.model.geo.synchronize()
687 gmsh.model.mesh.generate(3)
688
689 # Get mesh data for rotation
690 m = get_gmsh_entities()
691
692 # Rotate twice by 120 degrees to complete the triangle
693 for i in range(2):
694 angle = (i + 1) * 2*np.pi/3 # 120 degrees in radians
695 gmsh_transform_general(m, 1000*(i+1), 1000000*(i+1), 1000000*(i+1), angle)
696
697 # Remove duplicate nodes
698 gmsh.model.mesh.removeDuplicateNodes()
699
700 # Translate to specified center
701 gmsh_translate(xc)
702
703 else:
704 # For non-symmetric mesh, create full triangle directly
705 x_tri = get_ref_triangle_points(xc, r)
706
707 # Create points for the full triangle
708 points = []
709 points.append(gmsh.model.geo.addPoint(xc[0], xc[1], xc[2], h)) # Center point
710 for v in x_tri:
711 points.append(gmsh.model.geo.addPoint(v[0], v[1], v[2], h))
712
713 # Create lines connecting vertices
714 lines = []
715 for i in range(2):
716 lines.append(gmsh.model.geo.addLine(points[i+1], points[i+2]))
717 lines.append(gmsh.model.geo.addLine(points[3], points[1]))
718
719 # Create curve loop and surface
720 c1 = gmsh.model.geo.addCurveLoop(lines)
721 s1 = gmsh.model.geo.addPlaneSurface([c1])
722
723 # Synchronize
724 gmsh.model.geo.synchronize()
725
726 # Embed center point in surface
727 gmsh.model.mesh.embed(0, [points[0]], 2, s1)
728
729 # Generate mesh
730 gmsh.model.mesh.generate(3)
731
732 # Write output files
733 gmsh.write(filename + '.msh')
734 if vtk_out:
735 gmsh.write(filename + '.vtk')
736
737 gmsh.finalize()
738

References geom_util.get_gmsh_entities(), geom_util.get_ref_triangle_points(), geom_util.gmsh_transform_general(), and geom_util.gmsh_translate().

Here is the call graph for this function:

Variable Documentation

◆ a

gmsh_particles.a

Definition at line 1398 of file gmsh_particles.py.

◆ axis

gmsh_particles.axis

Definition at line 1446 of file gmsh_particles.py.

◆ bar_length

gmsh_particles.bar_length

Definition at line 1413 of file gmsh_particles.py.

◆ bar_width

gmsh_particles.bar_width

Definition at line 1412 of file gmsh_particles.py.

◆ center

gmsh_particles.center

Definition at line 1409 of file gmsh_particles.py.

◆ filename

gmsh_particles.filename

Definition at line 1363 of file gmsh_particles.py.

◆ h

gmsh_particles.h

Definition at line 1363 of file gmsh_particles.py.

◆ hole_Lx

gmsh_particles.hole_Lx

Definition at line 1440 of file gmsh_particles.py.

◆ hole_Ly

gmsh_particles.hole_Ly

Definition at line 1440 of file gmsh_particles.py.

◆ inner_radius

gmsh_particles.inner_radius

Definition at line 1411 of file gmsh_particles.py.

◆ inp_dir

str gmsh_particles.inp_dir = './'

Definition at line 1354 of file gmsh_particles.py.

◆ length

gmsh_particles.length

Definition at line 1446 of file gmsh_particles.py.

◆ Lx

gmsh_particles.Lx

Definition at line 1378 of file gmsh_particles.py.

◆ Ly

gmsh_particles.Ly

Definition at line 1378 of file gmsh_particles.py.

◆ Lz

gmsh_particles.Lz

Definition at line 1420 of file gmsh_particles.py.

◆ mesh

gmsh_particles.mesh

Definition at line 1444 of file gmsh_particles.py.

◆ mesh_size

gmsh_particles.mesh_size

Definition at line 1430 of file gmsh_particles.py.

◆ outer_radius

gmsh_particles.outer_radius

Definition at line 1410 of file gmsh_particles.py.

◆ points

gmsh_particles.points

Definition at line 1404 of file gmsh_particles.py.

◆ r

gmsh_particles.r

Definition at line 1363 of file gmsh_particles.py.

◆ R

gmsh_particles.R

Definition at line 1398 of file gmsh_particles.py.

◆ r_inner

gmsh_particles.r_inner

Definition at line 1435 of file gmsh_particles.py.

◆ r_outer

gmsh_particles.r_outer

Definition at line 1435 of file gmsh_particles.py.

◆ rx

gmsh_particles.rx

Definition at line 1368 of file gmsh_particles.py.

◆ ry

gmsh_particles.ry

Definition at line 1368 of file gmsh_particles.py.

◆ rz

gmsh_particles.rz

Definition at line 1425 of file gmsh_particles.py.

◆ symm_flag

int gmsh_particles.symm_flag = 0

Definition at line 1358 of file gmsh_particles.py.

◆ symmetric_mesh

gmsh_particles.symmetric_mesh

Definition at line 1363 of file gmsh_particles.py.

◆ test_meshes

list gmsh_particles.test_meshes = ['circle', 'ellipse', 'sphere', 'cuboid', 'ellipsoid', 'rectangle', 'hexagon', 'drum2d', 'triangle', 'polygon', 'cylindrical2d_wall', 'cylinder', 'annulus_circle', 'annulus_rectangle', 'open_rectangle', 'open_pipe']

Definition at line 1356 of file gmsh_particles.py.

◆ theta

gmsh_particles.theta = np.pi/6.

Definition at line 1397 of file gmsh_particles.py.

◆ True

gmsh_particles.True

Definition at line 1363 of file gmsh_particles.py.

◆ v1

list gmsh_particles.v1 = [0., 0., 0.]

Definition at line 1399 of file gmsh_particles.py.

◆ v2

list gmsh_particles.v2 = [R*np.cos(0.5*theta), -R*np.sin(0.5*theta), 0.]

Definition at line 1400 of file gmsh_particles.py.

◆ v3

list gmsh_particles.v3 = [R + a, 0., 0.]

Definition at line 1402 of file gmsh_particles.py.

◆ v4

list gmsh_particles.v4 = [R*np.cos(0.5*theta), R*np.sin(0.5*theta), 0.]

Definition at line 1401 of file gmsh_particles.py.

◆ vtk_out

gmsh_particles.vtk_out

Definition at line 1363 of file gmsh_particles.py.

◆ wall_thickness

gmsh_particles.wall_thickness

Definition at line 1447 of file gmsh_particles.py.

◆ width

gmsh_particles.width

Definition at line 1388 of file gmsh_particles.py.

◆ xc

gmsh_particles.xc

Definition at line 1363 of file gmsh_particles.py.