1
<?xml version="1.0" encoding="UTF-8" ?>
3
<name>p0p1cv with weak velBC navier stokes spatial convergence test using a rotated domain such as to have non cartesian aligned BC</name>
4
<owner userid="btollit"/>
6
<problem_definition length="short" nprocs="1">
8
fluidity -v1 -l MMS_A.flml
9
fluidity -v1 -l MMS_B.flml
10
fluidity -v1 -l MMS_C.flml
14
<variable name="ab_convergence_vel" language="python">
15
from fluidity_tools import stat_parser as stat
18
a_error_x = stat("MMS_A.stat")["NS"]["VectorAbsoluteDifference%1"]["l2norm"][-1]
19
b_error_x = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%1"]["l2norm"][-1]
20
a_error_y = stat("MMS_A.stat")["NS"]["VectorAbsoluteDifference%2"]["l2norm"][-1]
21
b_error_y = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%2"]["l2norm"][-1]
23
a_error_inf_x = stat("MMS_A.stat")["NS"]["VectorAbsoluteDifference%1"]["max"][-1]
24
b_error_inf_x = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%1"]["max"][-1]
25
a_error_inf_y = stat("MMS_A.stat")["NS"]["VectorAbsoluteDifference%2"]["max"][-1]
26
b_error_inf_y = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%2"]["max"][-1]
38
ab_ratio_x = a_error_x / b_error_x
39
ab_ratio_y = a_error_y / b_error_y
40
ab_ratio_inf_x = a_error_inf_x / b_error_inf_x
41
ab_ratio_inf_y = a_error_inf_y / b_error_inf_y
42
ab_convergence_vel = [[log(ab_ratio_x, 2), log(ab_ratio_inf_x, 2)], [log(ab_ratio_y, 2), log(ab_ratio_inf_y, 2)]]
44
<variable name="ab_convergence_p" language="python">
45
from fluidity_tools import stat_parser as stat
48
a_error = stat("MMS_A.stat")["NS"]["ScalarAbsoluteDifference"]["l2norm"][-1]
49
b_error = stat("MMS_B.stat")["NS"]["ScalarAbsoluteDifference"]["l2norm"][-1]
51
a_error_inf = stat("MMS_A.stat")["NS"]["ScalarAbsoluteDifference"]["max"][-1]
52
b_error_inf = stat("MMS_B.stat")["NS"]["ScalarAbsoluteDifference"]["max"][-1]
60
ab_ratio = a_error / b_error
61
ab_ratio_inf = a_error_inf / b_error_inf
62
ab_convergence_p = [log(ab_ratio, 2), log(ab_ratio_inf, 2)]
64
<variable name="bc_convergence_vel" language="python">
65
from fluidity_tools import stat_parser as stat
68
b_error_x = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%1"]["l2norm"][-1]
69
c_error_x = stat("MMS_C.stat")["NS"]["VectorAbsoluteDifference%1"]["l2norm"][-1]
70
b_error_y = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%2"]["l2norm"][-1]
71
c_error_y = stat("MMS_C.stat")["NS"]["VectorAbsoluteDifference%2"]["l2norm"][-1]
73
b_error_inf_x = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%1"]["max"][-1]
74
c_error_inf_x = stat("MMS_C.stat")["NS"]["VectorAbsoluteDifference%1"]["max"][-1]
75
b_error_inf_y = stat("MMS_B.stat")["NS"]["VectorAbsoluteDifference%2"]["max"][-1]
76
c_error_inf_y = stat("MMS_C.stat")["NS"]["VectorAbsoluteDifference%2"]["max"][-1]
88
bc_ratio_x = b_error_x / c_error_x
89
bc_ratio_y = b_error_y / c_error_y
90
bc_ratio_inf_x = b_error_inf_x / c_error_inf_x
91
bc_ratio_inf_y = b_error_inf_y / c_error_inf_y
92
bc_convergence_vel = [[log(bc_ratio_x, 2), log(bc_ratio_inf_x, 2)], [log(bc_ratio_y, 2), log(bc_ratio_inf_y, 2)]]
94
<variable name="bc_convergence_p" language="python">
95
from fluidity_tools import stat_parser as stat
98
b_error = stat("MMS_B.stat")["NS"]["ScalarAbsoluteDifference"]["l2norm"][-1]
99
c_error = stat("MMS_C.stat")["NS"]["ScalarAbsoluteDifference"]["l2norm"][-1]
101
b_error_inf = stat("MMS_B.stat")["NS"]["ScalarAbsoluteDifference"]["max"][-1]
102
c_error_inf = stat("MMS_C.stat")["NS"]["ScalarAbsoluteDifference"]["max"][-1]
110
bc_ratio = b_error / c_error
111
bc_ratio_inf = b_error_inf / c_error_inf
112
bc_convergence_p = [log(bc_ratio, 2), log(bc_ratio_inf, 2)]
114
<variable name="a_finish_time" language="python">
115
from fluidity_tools import stat_parser as stat
116
a_finish_time = stat("MMS_A.stat")["ElapsedTime"]["value"][-1]
118
<variable name="b_finish_time" language="python">
119
from fluidity_tools import stat_parser as stat
120
b_finish_time = stat("MMS_B.stat")["ElapsedTime"]["value"][-1]
122
<variable name="c_finish_time" language="python">
123
from fluidity_tools import stat_parser as stat
124
c_finish_time = stat("MMS_C.stat")["ElapsedTime"]["value"][-1]
126
<variable name="a_div" language="python">
127
from fluidity_tools import stat_parser as stat
128
a_div = max(abs(stat("MMS_A.stat")["NS"]["ControlVolumeDivergence"]["max"][-1]), abs(stat("MMS_A.stat")["NS"]["ControlVolumeDivergence"]["min"][-1]))
130
<variable name="b_div" language="python">
131
from fluidity_tools import stat_parser as stat
132
b_div = max(abs(stat("MMS_B.stat")["NS"]["ControlVolumeDivergence"]["max"][-1]), abs(stat("MMS_A.stat")["NS"]["ControlVolumeDivergence"]["min"][-1]))
134
<variable name="c_div" language="python">
135
from fluidity_tools import stat_parser as stat
136
c_div = max(abs(stat("MMS_C.stat")["NS"]["ControlVolumeDivergence"]["max"][-1]), abs(stat("MMS_A.stat")["NS"]["ControlVolumeDivergence"]["min"][-1]))
138
<variable name="a_final_change_vel" language="python">
140
from fluidity_tools import stat_parser as stat
141
vtu = vtktools.vtu("MMS_A_1.vtu")
142
dt = stat("MMS_A.stat")["dt"]["value"][-1]
143
a_final_change_vel = max(max(abs(vtu.GetVectorField("Velocity")[:,0]-vtu.GetVectorField("OldVelocity")[:,0])/dt), max(abs(vtu.GetVectorField("Velocity")[:,1]-vtu.GetVectorField("OldVelocity")[:,1])/dt))
145
<variable name="b_final_change_vel" language="python">
147
from fluidity_tools import stat_parser as stat
148
vtu = vtktools.vtu("MMS_B_1.vtu")
149
dt = stat("MMS_B.stat")["dt"]["value"][-1]
150
b_final_change_vel = max(max(abs(vtu.GetVectorField("Velocity")[:,0]-vtu.GetVectorField("OldVelocity")[:,0])/dt), max(abs(vtu.GetVectorField("Velocity")[:,1]-vtu.GetVectorField("OldVelocity")[:,1])/dt))
152
<variable name="c_final_change_vel" language="python">
154
from fluidity_tools import stat_parser as stat
155
vtu = vtktools.vtu("MMS_C_1.vtu")
156
dt = stat("MMS_C.stat")["dt"]["value"][-1]
157
c_final_change_vel = max(max(abs(vtu.GetVectorField("Velocity")[:,0]-vtu.GetVectorField("OldVelocity")[:,0])/dt), max(abs(vtu.GetVectorField("Velocity")[:,1]-vtu.GetVectorField("OldVelocity")[:,1])/dt))
159
<variable name="a_final_change_p" language="python">
161
from fluidity_tools import stat_parser as stat
162
vtu = vtktools.vtu("MMS_A_1.vtu")
163
dt = stat("MMS_A.stat")["dt"]["value"][-1]
164
a_final_change_p = max(abs(vtu.GetScalarField("Pressure")-vtu.GetScalarField("OldPressure"))/dt)
166
<variable name="b_final_change_p" language="python">
168
from fluidity_tools import stat_parser as stat
169
vtu = vtktools.vtu("MMS_B_1.vtu")
170
dt = stat("MMS_B.stat")["dt"]["value"][-1]
171
b_final_change_p = max(abs(vtu.GetScalarField("Pressure")-vtu.GetScalarField("OldPressure"))/dt)
173
<variable name="c_final_change_p" language="python">
175
from fluidity_tools import stat_parser as stat
176
vtu = vtktools.vtu("MMS_C_1.vtu")
177
dt = stat("MMS_C.stat")["dt"]["value"][-1]
178
c_final_change_p = max(abs(vtu.GetScalarField("Pressure")-vtu.GetScalarField("OldPressure"))/dt)
182
<test name="ab_convergence_vel_x: L2 order between 1.6 and 2.0" language="python">
183
assert(abs(ab_convergence_vel[0][0]-1.8) < 0.2)
185
<test name="ab_convergence_vel_y: L2 order between 1.6 and 2.0" language="python">
186
assert(abs(ab_convergence_vel[1][0]-1.8) < 0.2)
188
<test name="ab_convergence_p: L2 order between 2.0 and 2.2" language="python">
189
assert(abs(ab_convergence_p[0]-2.1) < 0.11)
191
<test name="bc_convergence_vel_x: L2 order between 0.8 and 1.0" language="python">
192
assert(abs(bc_convergence_vel[0][0]-0.9) < 0.1)
194
<test name="bc_convergence_vel_y: L2 order between 0.8 and 1.0" language="python">
195
assert(abs(bc_convergence_vel[1][0]-0.9) < 0.1)
197
<test name="bc_convergence_p: L2 order between 0.8 and 1.0" language="python">
198
assert(abs(bc_convergence_p[0]-0.9) < 0.1)
200
<test name="checking divergence in A vel with tolerance 1.0e-08" language="python">
201
assert(a_div < 1.E-8)
203
<test name="checking divergence in B vel with tolerance 1.0e-08" language="python">
204
assert(b_div < 1.E-8)
206
<test name="checking divergence in C vel with tolerance 1.0e-08" language="python">
207
assert(c_div < 1.E-8)
209
<test name="checking A finished in less than 50.0" language="python">
210
assert(a_finish_time-50.0 < 1.E-10)
212
<test name="checking B finished in less than 50.0" language="python">
213
assert(b_finish_time-50.0 < 1.E-10)
215
<test name="checking C finished in less than 50.0" language="python">
216
assert(c_finish_time-50.0 < 1.E-10)
218
<test name="checking steady state was reached in A vel with tolerance 1.0e-05" language="python">
219
assert(a_final_change_vel < 1.0e-05)
221
<test name="checking steady state was reached in B vel with tolerance 1.0e-05" language="python">
222
assert(b_final_change_vel < 1.0e-05)
224
<test name="checking steady state was reached in C vel with tolerance 1.0e-05" language="python">
225
assert(c_final_change_vel < 1.0e-05)
227
<test name="checking steady state was reached in A p with tolerance 1.0e-05" language="python">
228
assert(a_final_change_p < 1.0e-05)
230
<test name="checking steady state was reached in B p with tolerance 1.0e-05" language="python">
231
assert(b_final_change_p < 1.0e-05)
233
<test name="checking steady state was reached in C p with tolerance 1.0e-05" language="python">
234
assert(c_final_change_p < 1.0e-05)