-
Notifications
You must be signed in to change notification settings - Fork 1
/
mesh_generator.f90
185 lines (153 loc) · 9.62 KB
/
mesh_generator.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
!!----------------------------------------------------------------------------*
!! *** Program MESH_GENERATOR *** *
!! *
!! This module generates a basic square 2D traingular mesh. The *
!! geometry of this mesh is kept simple by placing nodes in a *
!! regular pattern, as demonstrated below (note the arrangment *
!! of the node and element IDs). *
!! *
!! 21._____22._____23._____24._____25. *
!! | / | / | / | / | *
!! | 29/ 25| 30/ 26| 31/ 27| 32/ 28| *
!! | / | / | / | / | *
!! 16._____17._____18._____19._____20. *
!! | / | / | / | / | *
!! | 21/ 17| 22/ 18| 23/ 19| 24/ 20| *
!! | / | / | / | / | *
!! 11._____12._____13._____14._____15. *
!! | / | / | / | / | *
!! | 13/ 9| 14/ 10| 15/ 11| 16/ 12| *
!! | / | / | / | / | *
!! 6.______7.______8.______9._____10. *
!! | / | / | / | / | *
!! | 5/ 1| 6/ 2| 7/ 3| 8/ 4| *
!! | / | / | / | / | *
!! 1.______2.______3.______4.______5. *
!! *
!! The size of the output square and its elements is determined *
!! by two inputs, box_size (the length of each side of the *
!! square) and edge_size (the length of a single elements *
!! boundary edge). Therefore, the above mesh would be generated *
!! by inputs that satisfy - box_size / edge_size = 4. *
!! *
!!----------------------------------------------------------------------------*
module mesh_generator
use, intrinsic :: iso_fortran_env
implicit none
contains
subroutine calculate_mesh_parameters(box_size, edge_size, num_edges_per_boundary, num_nodes, num_boundary_nodes, num_elements)
implicit none
integer(kind=int64), intent(in) :: box_size
real(kind=real64), intent(in) :: edge_size
integer(kind=int64), intent(out) :: num_edges_per_boundary, num_nodes, num_boundary_nodes, num_elements
num_edges_per_boundary = floor(box_size / edge_size)
num_nodes = (num_edges_per_boundary + 1)**2
num_boundary_nodes = (num_edges_per_boundary) * 4
num_elements = 2 * (num_edges_per_boundary)**2
! write(*,*) "** Mesh Parameters **"
! write(*,*) " num_edges_per_boundary:", num_edges_per_boundary
! write(*,*) " num_nodes:", num_nodes
! write(*,*) " num_boundary_nodes:", num_boundary_nodes
! write(*,*) " num_elements:", num_elements
end subroutine calculate_mesh_parameters
subroutine calculate_mesh(num_edges_per_boundary, num_nodes, num_elements, num_boundary_nodes, nodes, elements, boundary_edges)
implicit none
integer(kind=int64), intent(in) :: num_edges_per_boundary, num_nodes, num_boundary_nodes, num_elements
integer(kind=int64), dimension(3, num_elements), intent(inout) :: elements
integer(kind=int64), dimension(3, num_boundary_nodes), intent(inout) :: boundary_edges
real(kind=real64), dimension(2, num_nodes), intent(inout) :: nodes
integer :: num_nodes_per_boundary, bottom_left_node, counter, i, j
num_nodes_per_boundary = num_edges_per_boundary + 1
counter = 1
do i = 1, num_nodes_per_boundary
do j = 1, num_nodes_per_boundary
nodes(1, counter) = i ! x coordinate
nodes(2, counter) = j ! y coordinate
counter = counter + 1
end do
end do
counter = 1
do i = 1, num_edges_per_boundary
do j = 1, num_edges_per_boundary
bottom_left_node = j + (i - 1) * num_nodes_per_boundary
elements(1, counter) = bottom_left_node ! bottom left node
elements(2, counter) = bottom_left_node + 1 ! Next node anti-clockwise
elements(3, counter) = bottom_left_node + 1 + num_nodes_per_boundary ! Next node anti-clockwise
elements(1, counter + num_edges_per_boundary) = bottom_left_node
elements(2, counter + num_edges_per_boundary) = bottom_left_node + num_nodes_per_boundary + 1
elements(3, counter + num_edges_per_boundary) = bottom_left_node + num_nodes_per_boundary
counter = counter + 1
end do
counter = counter + num_edges_per_boundary
end do
! If we are along the bottom boundary
do i = 1, num_edges_per_boundary
! bottom boundary
boundary_edges(1, i) = i ! left node
boundary_edges(2, i) = i + 1 ! right node
boundary_edges(3, i) = i ! element
! right boundary
boundary_edges(1, i + num_edges_per_boundary) = i * num_nodes_per_boundary
boundary_edges(2, i + num_edges_per_boundary) = (i + 1) * num_nodes_per_boundary
boundary_edges(3, i + num_edges_per_boundary) = (2*i - 1) * num_edges_per_boundary
! top boundary
boundary_edges(1, i + num_edges_per_boundary * 2) = num_nodes - i + 1
boundary_edges(2, i + num_edges_per_boundary * 2) = num_nodes - i
boundary_edges(3, i + num_edges_per_boundary * 2) = num_elements - i + 1
! left boundary
boundary_edges(1, i + num_edges_per_boundary * 3) = (num_nodes_per_boundary - i) * num_nodes_per_boundary + 1
boundary_edges(2, i + num_edges_per_boundary * 3) = (num_nodes_per_boundary - 1 - i) * num_nodes_per_boundary + 1
boundary_edges(3, i + num_edges_per_boundary * 3) = num_elements - (num_edges_per_boundary - 1) - (2 * (i - 1) * num_edges_per_boundary)
end do
end subroutine calculate_mesh
subroutine write_mesh_to_file(num_nodes, num_elements, num_boundary_nodes, nodes, elements, boundary_edges)
implicit none
integer(kind=int64), intent(in) :: num_nodes, num_elements, num_boundary_nodes
integer(kind=int64), dimension(3, num_elements), intent(inout) :: elements
integer(kind=int64), dimension(3, num_boundary_nodes), intent(inout) :: boundary_edges
real(kind=real64), dimension(2, num_nodes), intent(inout) :: nodes
character*11 :: file_name
integer :: file_io
integer :: iostat
integer :: i, num_sets, num_dirichlet_boundary_conditions, num_neumann_boundary_conditions
! Baked in defaults
num_sets = 1
num_dirichlet_boundary_conditions = 1
num_neumann_boundary_conditions = 0
file_name = "square_mesh"
file_io = 100
! Write outpout
open (unit=file_io, &
file=file_name, &
status="replace", &
IOSTAT=iostat)
if( iostat .ne. 0) then
write(*,'(a)') ' *** Error when opening '//trim(file_name)
else
write(*,'(/,a)') ' *** '//trim(file_name)//' opened'
end if
write(file_io,*) "! num_nodes, num_elements, num_boundary_points, num_sets, num_dirichlet_boundary_conditions, num_neumann_boundary_conditions"
write(file_io,*) num_nodes, num_elements, num_boundary_nodes, num_sets, num_dirichlet_boundary_conditions, num_neumann_boundary_conditions
write(file_io,*) "! jb,vb(1,jb),vb(2,jb),vb(3,jb) - as many lines as num_sets"
write(file_io,*) 1, 1, 1, 1
write(file_io,*) "! jb,vb1(jb) - as many lines as num_dirichlet_boundary_conditions"
write(file_io,*) 1, 0
write(file_io,*) "! jb,vb2(jb) - as many lines as num_neumann_boundary_conditions"
write(file_io,*) "! jp,coordinates(1,jp),coordinates(2,jp) - as many lines as num_nodes"
do i = 1, num_nodes
write(file_io,*) i, nodes(1, i), nodes(2, i)
end do
write(file_io,*) "! je,element_to_node(1,je),element_to_node(2,je),element_to_node(3,je),vb_index(je) - as many lines as num_elements"
do i = 1, num_elements
write(file_io,*) i, elements(1, i), elements(2, i), elements(3, i), 1
end do
write(file_io,*) "! boundary_node_num(1,ib),boundary_node_num(2,ib) - as many lines as num_boundary_points"
do i = 1, num_boundary_nodes
write(file_io,*) i, 1
end do
write(file_io,*) "! num_side_nodes(1,ib),num_side_nodes(2,ib),num_side_nodes(3,ib),num_side_nodes(4,ib) - as many lines as num_boundary_points"
do i = 1, num_boundary_nodes
write(file_io,*) boundary_edges(1, i), boundary_edges(2, i), boundary_edges(3, i), 0
end do
end subroutine write_mesh_to_file
end module mesh_generator