forked from 2012ZGZYY/Dual_error_DG
-
Notifications
You must be signed in to change notification settings - Fork 0
/
indicator.cc
200 lines (176 loc) · 8 KB
/
indicator.cc
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
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/dofs/dof_handler.h>
#include "equation.h"
#include "claw.h"
using namespace dealii;
//-----------------------------------------------------------------------------
// compute KXRCF shock indicator based on density or energy
//-----------------------------------------------------------------------------
template <int dim>
void ConservationLaw<dim>::compute_shock_indicator ()
{
// If indicator type is "limiter" then mark all cells.
if(parameters.shock_indicator_type == Parameters::Limiter::limiter)
{
shock_indicator = 1e20;
}
else if(parameters.shock_indicator_type == Parameters::Limiter::u2)
{
compute_shock_indicator_u2();
}
else
{
compute_shock_indicator_kxrcf();
}
}
//-----------------------------------------------------------------------------
template <int dim>
void ConservationLaw<dim>::compute_shock_indicator_u2 ()
{
typename DoFHandler<dim>::active_cell_iterator
cell = dof_handler.begin_active(),
endc = dof_handler.end();
for(; cell != endc; ++cell)
{
bool u2 = test_u2(cell);
unsigned int c = cell_number(cell);
shock_indicator[c] = (u2) ? 0.0 : 1.0e20;
}
}
//-----------------------------------------------------------------------------
template <int dim>
void ConservationLaw<dim>::compute_shock_indicator_kxrcf ()
{
const unsigned int density_component = EulerEquations<dim>::density_component;
const unsigned int energy_component = EulerEquations<dim>::energy_component;
QGauss<dim-1> quadrature(fe.degree + 1);
FEFaceValues<dim> fe_face_values (mapping(), fe, quadrature,
update_values | update_normal_vectors);
FEFaceValues<dim> fe_face_values_nbr (mapping(), fe, quadrature,
update_values);
FESubfaceValues<dim> fe_subface_values (mapping(), fe, quadrature,
update_values | update_normal_vectors);
FESubfaceValues<dim> fe_subface_values_nbr (mapping(), fe, quadrature,
update_values);
unsigned int n_q_points = quadrature.size();
std::vector<double> face_values(n_q_points), face_values_nbr(n_q_points);
// select indicator variable
unsigned int component;
switch(parameters.shock_indicator_type)
{
case Parameters::Limiter::density:
component = density_component;
break;
case Parameters::Limiter::energy:
component = energy_component;
break;
default:
component = 0;
AssertThrow(false, ExcNotImplemented());
}
const FEValuesExtractors::Scalar variable (component);
typename DoFHandler<dim>::active_cell_iterator
cell = dof_handler.begin_active(),
endc = dof_handler.end();
jump_ind_min = 1.0e20;
jump_ind_max = 0.0;
jump_ind_avg = 0.0;
for(; cell != endc; ++cell)
{
unsigned int c = cell_number(cell);
double& cell_shock_ind = shock_indicator (c);
double& cell_jump_ind = jump_indicator (c);
cell_shock_ind = 0;
cell_jump_ind = 0;
double inflow_measure = 0;
// velocity based on cell average. we use this to determine inflow/outflow
// parts of the cell boundary.
Point<dim> vel;
for(unsigned int i=0; i<dim; ++i)
vel(i) = cell_average[c][i] / cell_average[c][density_component];
for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
if (cell->at_boundary(f) == false)
{
if ((cell->neighbor(f)->level() == cell->level()) &&
(cell->neighbor(f)->has_children() == false))
{
fe_face_values.reinit(cell, f);
fe_face_values_nbr.reinit(cell->neighbor(f), cell->neighbor_of_neighbor(f));
fe_face_values[variable].get_function_values(current_solution, face_values);
fe_face_values_nbr[variable].get_function_values(current_solution, face_values_nbr);
for(unsigned int q=0; q<n_q_points; ++q)
{
int inflow_status = (vel * fe_face_values.normal_vector(q) < 0);
cell_shock_ind += inflow_status *
(face_values[q] - face_values_nbr[q]) *
fe_face_values.JxW(q);
cell_jump_ind += std::pow(face_values[q] - face_values_nbr[q], 2) *
fe_face_values.JxW(q);
inflow_measure += inflow_status * fe_face_values.JxW(q);
}
}
else if ((cell->neighbor(f)->level() == cell->level()) &&
(cell->neighbor(f)->has_children() == true))
{
for (unsigned int subface=0; subface<cell->face(f)->n_children(); ++subface)
{
fe_subface_values.reinit (cell, f, subface);
fe_face_values_nbr.reinit (cell->neighbor_child_on_subface (f, subface),
cell->neighbor_of_neighbor(f));
fe_subface_values[variable].get_function_values(current_solution, face_values);
fe_face_values_nbr[variable].get_function_values(current_solution, face_values_nbr);
for(unsigned int q=0; q<n_q_points; ++q)
{
int inflow_status = (vel * fe_subface_values.normal_vector(q) < 0);
cell_shock_ind += inflow_status *
(face_values[q] - face_values_nbr[q]) *
fe_subface_values.JxW(q);
cell_jump_ind += std::pow(face_values[q] - face_values_nbr[q], 2) *
fe_face_values.JxW(q);
inflow_measure += inflow_status * fe_subface_values.JxW(q);
}
}
}
else if (cell->neighbor_is_coarser(f))
{
fe_face_values.reinit(cell, f);
fe_subface_values_nbr.reinit (cell->neighbor(f),
cell->neighbor_of_coarser_neighbor(f).first,
cell->neighbor_of_coarser_neighbor(f).second);
fe_face_values[variable].get_function_values(current_solution, face_values);
fe_subface_values_nbr[variable].get_function_values(current_solution, face_values_nbr);
for(unsigned int q=0; q<n_q_points; ++q)
{
int inflow_status = (vel * fe_face_values.normal_vector(q) < 0);
cell_shock_ind += inflow_status *
(face_values[q] - face_values_nbr[q]) *
fe_face_values.JxW(q);
cell_jump_ind += std::pow(face_values[q] - face_values_nbr[q], 2) *
fe_face_values.JxW(q);
inflow_measure += inflow_status * fe_face_values.JxW(q);
}
}
}
else
{
// Boundary face
// We dont do anything here since we assume solution is constant near
// boundary.
}
// normalized shock indicator
double cell_norm = cell_average[c][component];
double denominator = std::pow(cell->diameter(), 0.5*(fe.degree+1)) *
inflow_measure *
cell_norm;
cell_shock_ind = std::fabs(cell_shock_ind) / denominator;
double dx = cell->diameter() / std::sqrt(1.0*dim);
// TODO: normalization needs to be done properly
cell_jump_ind = std::sqrt( cell_jump_ind / (4.0*dx) ) * cell->diameter();
jump_ind_min = std::min(jump_ind_min, cell_jump_ind);
jump_ind_max = std::max(jump_ind_max, cell_jump_ind);
jump_ind_avg += cell_jump_ind;
}
jump_ind_avg /= triangulation.n_active_cells();
}
template class ConservationLaw<2>;