-
Notifications
You must be signed in to change notification settings - Fork 95
/
Copy pathtest_grid_dofhandler_vtk.jl
166 lines (133 loc) · 5.51 KB
/
test_grid_dofhandler_vtk.jl
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
# to test vtk-files
OVERWRITE_CHECKSUMS = false
checksums_file = joinpath(dirname(@__FILE__), "checksums.sha1")
checksum_list = read(checksums_file, String)
if OVERWRITE_CHECKSUMS
csio = open(checksums_file, "w")
else
csio = open(checksums_file, "r")
end
@testset "Grid, DofHandler, vtk" begin
for (celltype, dim) in ((Line, 1),
(QuadraticLine, 1),
(Quadrilateral, 2),
(QuadraticQuadrilateral, 2),
(Triangle, 2),
(QuadraticTriangle, 2),
(Hexahedron, 3),
(Tetrahedron, 3))
# create test grid, do some operations on it and then test
# the resulting sha1 of the stored vtk file
# after manually checking the exported vtk
nels = ntuple(x->5, dim)
right = Vec{dim, Float64}(ntuple(x->1.5, dim))
left = -right
grid = generate_grid(celltype, nels, left, right)
transform!(grid, x-> 2x)
radius = 2*1.5
addcellset!(grid, "cell-1", [1,])
addcellset!(grid, "middle-cells", x -> norm(x) < radius)
addfaceset!(grid, "middle-faceset", x -> norm(x) < radius)
addfaceset!(grid, "right-faceset", getfaceset(grid, "right"))
addnodeset!(grid, "middle-nodes", x -> norm(x) < radius)
gridfilename = "grid-$(Ferrite.celltypes[celltype])"
vtk_grid(gridfilename, grid) do vtk
vtk_cellset(vtk, grid, "cell-1")
vtk_cellset(vtk, grid, "middle-cells")
vtk_nodeset(vtk, grid, "middle-nodes")
end
# test the sha of the file
sha = bytes2hex(open(SHA.sha1, gridfilename*".vtu"))
if OVERWRITE_CHECKSUMS
write(csio, sha, "\n")
else
@test chomp(readline(csio)) == sha
rm(gridfilename*".vtu")
end
# Create a DofHandler, add some things, write to file and
# then check the resulting sha
dofhandler = DofHandler(grid)
push!(dofhandler, :temperature, 1)
push!(dofhandler, :displacement, 3)
close!(dofhandler)
ch = ConstraintHandler(dofhandler)
dbc = Dirichlet(:temperature, union(getfaceset(grid, "left"), getfaceset(grid, "right-faceset")), (x,t)->1)
add!(ch, dbc)
dbc = Dirichlet(:temperature, getfaceset(grid, "middle-faceset"), (x,t)->4)
add!(ch, dbc)
for d in 1:dim
dbc = Dirichlet(:displacement, union(getfaceset(grid, "left")), (x,t) -> d, d)
add!(ch, dbc)
end
close!(ch)
update!(ch, 0.0)
Random.seed!(1234)
u = rand(ndofs(dofhandler))
apply!(u, ch)
dofhandlerfilename = "dofhandler-$(Ferrite.celltypes[celltype])"
vtk_grid(dofhandlerfilename, dofhandler) do vtk
vtk_point_data(vtk, ch)
vtk_point_data(vtk, dofhandler, u)
end
# test the sha of the file
sha = bytes2hex(open(SHA.sha1, dofhandlerfilename*".vtu"))
if OVERWRITE_CHECKSUMS
write(csio, sha, "\n")
else
@test chomp(readline(csio)) == sha
rm(dofhandlerfilename*".vtu")
end
end
end # of testset
if OVERWRITE_CHECKSUMS
close(csio)
end
# right = Vec{2, Float64}(ntuple(x->1.5, dim))
# left = -right
@testset "Grid utils" begin
grid = Ferrite.generate_grid(QuadraticQuadrilateral, (1, 1), Vec((0.,0.)), Vec((1.,1.)))
addcellset!(grid, "cell_set", [1]);
node_set = Set(1:getnnodes(grid))
addnodeset!(grid, "node_set", node_set)
@test getnodesets(grid) == Dict("node_set" => node_set)
@test getnodes(grid, [1]) == [getnodes(grid, 1)] # untested
@test length(getnodes(grid, "node_set")) == 9
@test collect(getcoordinates(getnodes(grid, 5)).data) ≈ [0.5, 0.5]
@test getcells(grid, "cell_set") == [getcells(grid, 1)]
f(x) = Tensor{1,1,Float64}((1 + x[1]^2 + 2x[2]^2, ))
values = compute_vertex_values(grid, f)
@test f([0.0, 0.0]) == values[1]
@test f([0.5, 0.5]) == values[5]
@test f([1.0, 1.0]) == values[9]
@test compute_vertex_values(grid, collect(1:9), f) == values
# Can we test this in a better way? The set makes the order random.
@test length(compute_vertex_values(grid, "node_set", f)) == 9
# CellIterator on a grid without DofHandler
grid = generate_grid(Triangle, (4,4))
n = 0
ci = CellIterator(grid)
@test length(ci) == getncells(grid)
for c in ci
getcoordinates(c)
getnodes(c)
n += cellid(c)
end
@test n == div(getncells(grid)*(getncells(grid) + 1), 2)
end
@testset "Grid sets" begin
grid = Ferrite.generate_grid(Hexahedron, (1, 1, 1), Vec((0.,0., 0.)), Vec((1.,1.,1.)))
#Test manual add
addcellset!(grid, "cell_set", [1]);
addnodeset!(grid, "node_set", [1])
addfaceset!(grid, "face_set", [FaceIndex(1,1)])
addedgeset!(grid, "edge_set", [EdgeIndex(1,1)])
addvertexset!(grid, "vert_set", [VertexIndex(1,1)])
#Test function add
addfaceset!(grid, "left_face", (x)-> x[1] ≈ 0.0)
addedgeset!(grid, "left_lower_edge", (x)-> x[1] ≈ 0.0 && x[3] ≈ 0.0)
addvertexset!(grid, "left_corner", (x)-> x[1] ≈ 0.0 && x[2] ≈ 0.0 && x[3] ≈ 0.0)
@test 1 in Ferrite.getnodeset(grid, "node_set")
@test FaceIndex(1,5) in getfaceset(grid, "left_face")
@test EdgeIndex(1,4) in getedgeset(grid, "left_lower_edge")
@test VertexIndex(1,1) in getvertexset(grid, "left_corner")
end