forked from uw-comphys/opencmp
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy patherror_analysis.py
executable file
·150 lines (120 loc) · 6.45 KB
/
error_analysis.py
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
"""
Copyright 2021 the authors (see AUTHORS file for full list)
This file is part of OpenCMP.
OpenCMP is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation, either version 2.1 of the License, or
(at your option) any later version.
OpenCMP is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with OpenCMP. If not, see <https://www.gnu.org/licenses/>.
"""
from config_functions import ConfigParser
from helpers.error import norm
from solvers import Solver
from ngsolve import GridFunction
import tabulate
import math
def h_convergence(config: ConfigParser, solver: Solver, sol: GridFunction, var: str) -> None:
"""
Function to check h (mesh element size) convergence and print results.
Args:
config: Config file from which to grab
solver: The solver used
sol: Gridfunction that contains the current solution
var: The variable of interest.
"""
num_refinements = config.get_item(['ERROR ANALYSIS', 'num_refinements'], int)
average_lst = config.get_list(['ERROR ANALYSIS', 'error_average'], str, quiet=True)
component = solver.model.model_components[var]
average = component in average_lst
# Reload the model's mesh and finite element space so convergence tests can be chained and don't affect each other.
solver.model.load_mesh_fes(mesh=True, fes=True)
# First solve used the default settings.
if component is None:
err = norm('l2_norm', sol, solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes, average)
else:
err = norm('l2_norm', sol.components[component], solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes.components[component], average)
# Track the convergence information.
num_dofs_lst = [solver.model.fes.ndof]
error_lst = [err]
# Then run through a series of mesh refinements and resolve on each
# refined mesh.
for n in range(num_refinements):
solver.model.mesh.Refine()
solver.model.fes.Update()
solver.reset_model()
sol = solver.solve()
if component is None:
err = norm('l2_norm', sol, solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes, average)
else:
err = norm('l2_norm', sol.components[component], solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes.components[component], average)
num_dofs_lst.append(solver.model.fes.ndof)
error_lst.append(err)
print('L2 norm at refinement {0}: {1}'.format(n, err))
# Display the results nicely.
convergence_table = [['Refinement Level', 'DOFs', 'Error', 'Convergence Rate']]
convergence_table.append([1, num_dofs_lst[0], error_lst[0], 0])
for n in range(num_refinements):
ref_level = '1/{}'.format(int(2 ** (n + 1)))
convergence_rate = math.log(error_lst[n] / error_lst[n + 1]) / \
(math.log(num_dofs_lst[n + 1] / (num_dofs_lst[n] * 2.0 ** (solver.model.mesh.dim - 1))))
convergence_table.append([ref_level, num_dofs_lst[n + 1], error_lst[n + 1], convergence_rate])
print(tabulate.tabulate(convergence_table, headers='firstrow', floatfmt=('.1f', '.1f', '.3e', '.2f')))
def p_convergence(config: ConfigParser, solver: Solver, sol: GridFunction, var: str) -> None:
"""
Function to check p (interpolat polynomial order) convergence and print results.
Args:
config: Config file from which to grab
solver: The solver used
sol: Gridfunction that contains the current solution
var: The variable of interest.
"""
num_refinements = config.get_item(['ERROR ANALYSIS', 'num_refinements'], int)
average_lst = config.get_list(['ERROR ANALYSIS', 'error_average'], str, quiet=True)
component = solver.model.model_components[var]
average = component in average_lst
# Reload the model's mesh and finite element space so convergence tests can be chained and don't affect each other.
solver.model.load_mesh_fes(mesh=True, fes=True)
# First solve used the default settings.
if component is None:
err = norm('l2_norm', sol, solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes, average)
else:
err = norm('l2_norm', sol.components[component], solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes.components[component], average)
# Track the convergence information.
num_dofs_lst = [solver.model.fes.ndof]
interp_ord_lst = [solver.model.interp_ord[var]]
error_lst = [err]
# Then run through a series of interpolant refinements.
for n in range(num_refinements):
solver.model.interp_ord = {key: val + 1 for key, val in solver.model.interp_ord.items()}
solver.model.load_mesh_fes(mesh=False, fes=True)
solver.reset_model()
sol = solver.solve()
if component is None:
err = norm('l2_norm', sol, solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes, average)
else:
err = norm('l2_norm', sol.components[component], solver.model.ref_sol['ref_sols'][var],
solver.model.mesh, solver.model.fes.components[component], average)
num_dofs_lst.append(solver.model.fes.ndof)
interp_ord_lst.append(solver.model.interp_ord[var])
error_lst.append(err)
print('L2 norm at refinement {0}: {1}'.format(n, err))
# Display the results nicely.
convergence_table = [['Interpolant Order', 'DOFs', 'Error', 'Convergence Rate']]
convergence_table.append([interp_ord_lst[0], num_dofs_lst[0], error_lst[0], 0])
for n in range(num_refinements):
# TODO: Not sure if this is the correct equation for p-convergence.
convergence_rate = math.log(error_lst[n] / error_lst[n + 1]) / math.log(num_dofs_lst[n + 1] / num_dofs_lst[n])
convergence_table.append([interp_ord_lst[n + 1], num_dofs_lst[n + 1], error_lst[n + 1], convergence_rate])
print(tabulate.tabulate(convergence_table, headers='firstrow', floatfmt=['.1f', '.1f', '.3e', '.2f']))