2 #include "l2_proj_based_selector.h"
8 namespace RefinementSelectors
10 template<
typename Scalar>
14 if (user_shapeset !=
nullptr)
16 this->warn(
"Warning: The user shapeset provided for the selector has to have a correct copy constructor implemented.");
17 this->warn(
"Warning: The functionality for cloning user shapeset is to be implemented yet.");
21 template<
typename Scalar>
24 if (!this->user_shapeset)
25 delete this->shapeset;
28 template<
typename Scalar>
31 int max_element_order = (20 - element->
iro_cache) / 2 - 1;
33 max_order_ = max_element_order;
35 max_order_ = std::min(this->max_order, max_element_order);
39 template<
typename Scalar>
45 while (!done && inx_trf < H2D_TRF_NUM)
48 const Trf& trf = trfs[inx_trf];
49 std::vector<typename ProjBasedSelector<Scalar>::TrfShapeExp>& trf_svals = svals[inx_trf];
52 trf_svals.resize(max_shape_inx + 1);
55 const int num_shapes = (int)shapes.size();
56 for (
int i = 0; i < num_shapes; i++)
58 int inx_shape = shapes[i].inx;
62 shape_exp.allocate(H2D_L2FE_NUM, num_gip_points);
65 for (
int k = 0; k < num_gip_points; k++)
68 double ref_x = gip_points[k][
H2D_GIP2D_X] * trf.m[0] + trf.
t[0];
69 double ref_y = gip_points[k][
H2D_GIP2D_Y] * trf.m[1] + trf.
t[1];
72 shape_exp[H2D_L2FE_VALUE][k] = this->shapeset->get_fn_value(inx_shape, ref_x, ref_y, 0, mode);
77 if (inx_trf == H2D_TRF_IDENTITY)
82 if (inx_trf >= num_noni_trfs)
83 inx_trf = H2D_TRF_IDENTITY;
87 throw Exceptions::Exception(
"All transformation processed but identity transformation not found.");
90 template<
typename Scalar>
94 precalc_shapes(gip_points, num_gip_points, trfs, num_noni_trfs, shapes, max_shape_inx, svals, mode);
97 const int num_shapes = (int)shapes.size();
98 for (
int i = 0; i < num_shapes; i++)
100 const int inx_shape_i = shapes[i].inx;
103 for (
int j = 0; j < i; j++)
105 const int inx_shape_j = shapes[j].inx;
108 double product = 0.0;
109 for (
int k = 0; k < num_gip_points; k++)
112 sum += svals[H2D_TRF_IDENTITY][inx_shape_i][H2D_L2FE_VALUE][k] * svals[H2D_TRF_IDENTITY][inx_shape_j][H2D_L2FE_VALUE][k];
119 while (!done && inx_trf < H2D_TRF_NUM)
122 for (
int k = 0; k < num_gip_points; k++)
124 svals[inx_trf][inx_shape_i][H2D_L2FE_VALUE][k] -= product * svals[inx_trf][inx_shape_j][H2D_L2FE_VALUE][k];
128 if (inx_trf == H2D_TRF_IDENTITY)
133 if (inx_trf >= num_noni_trfs)
134 inx_trf = H2D_TRF_IDENTITY;
139 throw Exceptions::Exception(
"All transformation processed but identity transformation not found.");
144 double norm_squared = 0.0;
145 for (
int k = 0; k < num_gip_points; k++)
148 sum += sqr(svals[H2D_TRF_IDENTITY][inx_shape_i][H2D_L2FE_VALUE][k]);
151 double norm = sqrt(norm_squared);
152 if (!finite(1 / norm))
153 throw Exceptions::Exception(
"Norm (%g) is almost zero.", norm);
158 while (!done && inx_trf < H2D_TRF_NUM)
161 for (
int k = 0; k < num_gip_points; k++)
163 svals[inx_trf][inx_shape_i][H2D_L2FE_VALUE][k] /= norm;
167 if (inx_trf == H2D_TRF_IDENTITY)
172 if (inx_trf >= num_noni_trfs)
173 inx_trf = H2D_TRF_IDENTITY;
178 throw Exceptions::Exception(
"All transformation processed but identity transformation not found.");
182 template<
typename Scalar>
186 rval[inx_son][H2D_L2FE_VALUE] = malloc_with_check<Scalar>(num_gip);
187 memcpy(rval[inx_son][H2D_L2FE_VALUE], rsln->
get_fn_values(0), num_gip *
sizeof(Scalar));
190 template<
typename Scalar>
193 free_with_check(rval[inx_son][H2D_L2FE_VALUE]);
196 template<
typename Scalar>
198 const int* shape_inx,
const int num_shapes,
ElementMode2D mode)
201 double** matrix = new_matrix<double>(num_shapes, num_shapes);
205 for (
int i = 0; i < num_shapes; i++, inx_row += num_shapes)
207 double* matrix_row = matrix[i];
208 int shape0_inx = shape_inx[i];
209 for (
int k = 0; k < num_shapes; k++)
211 int shape1_inx = shape_inx[k];
214 for (
int j = 0; j < num_gip_points; j++)
217 double value0 = this->shapeset->get_value(
H2D_FEI_VALUE, shape0_inx, gip_x, gip_y, 0, mode);
218 double value1 = this->shapeset->get_value(
H2D_FEI_VALUE, shape1_inx, gip_x, gip_y, 0, mode);
220 value += gip_points[j][
H2D_GIP2D_W] * (value0*value1);
223 matrix_row[k] = value;
230 template<
typename Scalar>
233 Scalar total_value = 0;
234 for (
int gip_inx = 0; gip_inx < sub_gip.
num_gip_points; gip_inx++)
237 double3 &gip_pt = sub_gip.
gip_points[gip_inx];
240 double shape_value = sub_shape.
svals[H2D_L2FE_VALUE][gip_inx];
244 ref_value = rval[son][H2D_L2FE_VALUE][gip_inx];
247 Scalar value = (shape_value * ref_value);
254 template<
typename Scalar>
257 double total_error_squared = 0;
258 for (
int gip_inx = 0; gip_inx < sub_gip.
num_gip_points; gip_inx++)
261 double3 &gip_pt = sub_gip.
gip_points[gip_inx];
264 Scalar proj_value = 0;
265 for (
int i = 0; i < elem_proj.
num_shapes; i++)
268 proj_value += elem_proj.
shape_coeffs[i] * elem_proj.
svals[shape_inx][H2D_L2FE_VALUE][gip_inx];
272 Scalar ref_value = rval[son][H2D_L2FE_VALUE][gip_inx];
275 double error_squared = sqr(proj_value - ref_value);
277 total_error_squared += gip_pt[
H2D_GIP2D_W] * error_squared;
280 return total_error_squared;
virtual void get_current_order_range(Element *element, int &min_order, int &max_order)
Sets OptimumSelector::current_max_order and OptimumSelector::current_min_order.
virtual void precalc_shapes(const double3 *gip_points, const int num_gip_points, const Trf *trfs, const int num_noni_trfs, const std::vector< typename OptimumSelector< Scalar >::ShapeInx > &shapes, const int max_shape_inx, typename ProjBasedSelector< Scalar >::TrfShape &svals, ElementMode2D mode)
Calculates values of shape function at GIP for all transformations.
A transform shaped function expansions.
Stores one element of a mesh.
int * shape_inxs
Used shape indices.
#define H2DRS_MAX_ORDER
A maximum order suported by refinement selectors.
virtual void precalc_ref_solution(int inx_son, MeshFunction< Scalar > *rsln, Element *element, int intr_gip_order, Scalar *rval[H2D_MAX_ELEMENT_SONS][MAX_NUMBER_FUNCTION_VALUES_FOR_SELECTORS])
Returns an array of values of the reference solution at integration points.
A shape function on subdomain of an element.
Represents a function defined on a mesh.
Quad2D * get_quad_2d() const
Returns the current quadrature points.
Integration points in the reference domain of an element of a candidate.
Projection of an element of a candidate.
L2 shapeset - products of legendre polynomials.
double3 * gip_points
Integration points and weights. The first index is an index of an integration point, the second index is defined through the enum GIP2DIndices.
A general projection-based selector.
int num_gip_points
A number of integration points.
int num_shapes
A number of used shape indices.
CandList
Predefined list of candidates.
Index of a function value f.
unsigned short iro_cache
Increase in integration order, see RefMap::calc_inv_ref_order()
virtual void free_ref_solution_data(int inx_son, Scalar *rval[H2D_MAX_ELEMENT_SONS][MAX_NUMBER_FUNCTION_VALUES_FOR_SELECTORS])
Frees the data allocated in precalc_ref_solution.
Scalar * shape_coeffs
Coefficients of shape indices of a projection.
double2 t
The 2x2 diagonal transformation matrix.
std::vector< TrfShapeExp > & svals
A precalculated shape-function values. Empty is not defined.
TrfShapeExp & svals
Evaluate values of a shape function. If TrfShapeExp::empty(), no precalculated values are available...
A transformation from a reference domain of a subdomain to a reference domain of an element of a cand...
L2ProjBasedSelector(CandList cand_list=H2D_HP_ANISO, int max_order=H2DRS_DEFAULT_ORDER, L2Shapeset *user_shapeset=nullptr)
Constructor.
#define H2D_MAX_ELEMENT_SONS
Macros.
virtual double evaluate_error_squared_subdomain(Element *sub_elem, const typename ProjBasedSelector< Scalar >::ElemGIP &sub_gip, int son, const typename ProjBasedSelector< Scalar >::ElemSubTrf &sub_trf, const typename ProjBasedSelector< Scalar >::ElemProj &elem_proj, Scalar *rval[H2D_MAX_ELEMENT_SONS][MAX_NUMBER_FUNCTION_VALUES_FOR_SELECTORS])
Evaluates an squared error of a projection of an element of a candidate onto subdomains.
virtual const Scalar * get_fn_values(int component=0) const
Returns function values.
~L2ProjBasedSelector()
Destructor.
virtual void precalc_ortho_shapes(const double3 *gip_points, const int num_gip_points, const Trf *trfs, const int num_noni_trfs, const std::vector< typename OptimumSelector< Scalar >::ShapeInx > &shapes, const int max_shape_inx, typename ProjBasedSelector< Scalar >::TrfShape &svals, ElementMode2D mode)
Calculates values of orthogonalized shape function at GIP for all transformations.
A projection-based selector for L2 space.
#define H2DRS_DEFAULT_ORDER
A default order. Used to indicate an unkonwn order or a maximum support order.
virtual Scalar evaluate_rhs_subdomain(Element *sub_elem, const typename ProjBasedSelector< Scalar >::ElemGIP &sub_gip, int son, const typename ProjBasedSelector< Scalar >::ElemSubTrf &sub_trf, const typename ProjBasedSelector< Scalar >::ElemSubShapeFunc &sub_shape, Scalar *rval[H2D_MAX_ELEMENT_SONS][MAX_NUMBER_FUNCTION_VALUES_FOR_SELECTORS])
Evaluates a value of the right-hande side in a subdomain.
virtual double ** build_projection_matrix(double3 *gip_points, int num_gip_points, const int *shape_inx, const int num_shapes, ElementMode2D mode)
Builds projection matrix using a given set of shapes.