92 #ifndef GETFEM_INTEGRATION_H__
93 #define GETFEM_INTEGRATION_H__
116 mutable std::vector<long_scalar_type> int_monomials;
117 mutable std::vector< std::vector<long_scalar_type> > int_face_monomials;
122 dim_type
dim(
void)
const {
return cvs->dim(); }
132 long_scalar_type
int_poly(
const base_poly &P)
const;
143 typedef std::shared_ptr<const poly_integration> ppoly_integration;
151 class integration_method;
152 typedef std::shared_ptr<const integration_method> pintegration_method;
155 class approx_integration {
159 bgeot::pconvex_ref cvr;
160 bgeot::pstored_point_tab pint_points;
161 std::vector<scalar_type> int_coeffs;
162 std::vector<size_type> repartition;
165 std::vector<PT_TAB> pt_to_store;
167 bool built_on_the_fly;
172 dim_type dim(
void)
const {
return cvr->structure()->dim(); }
173 size_type nb_points(
void)
const {
return int_coeffs.size(); }
175 size_type nb_points_on_convex(
void)
const {
return repartition[0]; }
178 {
return repartition[f+1] - repartition[f]; }
180 {
return repartition[f]; }
183 {
return cvr->structure()->nb_faces(); }
184 bool is_built_on_the_fly(
void)
const {
return built_on_the_fly; }
185 void set_built_on_the_fly(
void) { built_on_the_fly =
true; }
189 bgeot::pconvex_ref ref_convex(
void)
const {
return cvr; }
191 const std::vector<size_type> &repart(
void)
const {
return repartition; }
196 bgeot::pstored_point_tab pintegration_points(
void)
const
197 {
return pint_points; }
199 const base_node &point(
size_type i)
const
200 {
return (*pint_points)[i]; }
204 {
return (*pint_points)[repartition[f] + i]; }
206 const std::vector<scalar_type> &integration_coefficients(
void)
const
207 {
return int_coeffs; }
209 scalar_type coeff(
size_type i)
const {
return int_coeffs[i]; }
212 {
return int_coeffs[repartition[f] + i]; }
214 void add_point(
const base_node &pt, scalar_type w,
216 bool include_empty =
false);
217 void add_point_norepeat(
const base_node &pt, scalar_type w,
219 void add_point_full_symmetric(base_node pt, scalar_type w);
220 void add_method_on_face(pintegration_method ppi,
short_type f);
221 void valid_method(
void);
223 approx_integration(
void) : valid(false), built_on_the_fly(false) { }
224 approx_integration(bgeot::pconvex_ref cr)
225 : cvr(cr), repartition(cr->structure()->nb_faces()+1),
226 pt_to_store(cr->structure()->nb_faces()+1), valid(false),
227 built_on_the_fly(false)
228 { std::fill(repartition.begin(), repartition.end(), 0); }
229 virtual ~approx_integration() {}
232 typedef std::shared_ptr<const approx_integration> papprox_integration;
244 ppoly_integration ppi;
245 papprox_integration pai;
247 void remove() { pai = papprox_integration(); ppi = ppoly_integration(); }
251 const papprox_integration &approx_method(
void)
const {
return pai; }
252 const ppoly_integration &exact_method(
void)
const {
return ppi; }
254 void set_approx_method(papprox_integration paii)
255 { remove(); pai = paii; im_type = IM_APPROX; }
256 void set_exact_method(ppoly_integration ppii)
257 { remove(); ppi = ppii; im_type = IM_EXACT; }
259 bgeot::pstored_point_tab pintegration_points(
void)
const {
260 if (type() == IM_EXACT) {
262 std::vector<base_node> spt(1); spt[0] = base_node(n);
263 return store_point_tab(spt);
265 else if (type() == IM_APPROX)
266 return pai->pintegration_points();
267 else GMM_ASSERT1(
false,
"IM_NONE has no points");
275 case IM_EXACT:
return ppi->structure();
276 case IM_APPROX:
return pai->structure();
277 case IM_NONE: GMM_ASSERT1(
false,
"IM_NONE has no structure");
278 default: GMM_ASSERT3(
false,
"");
284 DAL_STORED_OBJECT_DEBUG_CREATED(
this,
"Exact integration method");
285 ppi = p; im_type = IM_EXACT;
289 DAL_STORED_OBJECT_DEBUG_CREATED(
this,
"Approximate integration method");
290 pai = p; im_type = IM_APPROX;
294 DAL_STORED_OBJECT_DEBUG_CREATED(
this,
"Integration method");
299 { DAL_STORED_OBJECT_DEBUG_DESTROYED(
this,
"Integration method"); }
311 bool throw_if_not_found =
true);
339 pintegration_method
im_none(
void);
343 papprox_integration composite_approx_int_method(
const bgeot::mesh_precomposite &mp,
345 bgeot::pconvex_ref cr);
349 scalar_type test_integration_error(papprox_integration pim, dim_type order);
351 papprox_integration get_approx_im_or_fail(pintegration_method pim);
356 typedef dal::naming_system<integration_method>::param_list im_param_list;
358 void add_integration_name(std::string name,
Geometric transformations on convexes.
Structure which dynamically collects points identifying points that are nearer than a certain very sm...
Handle composite polynomials.
Store a set of points, identifying points that are nearer than a certain very small distance.
Vector of integer (16 bits type) which represent the powers of a monomial.
Associate a name to a method descriptor and store method descriptors.
base class for static stored objects
this structure is not intended to be used directly.
Describe an integration method linked to a mesh.
Description of an exact integration of polynomials.
long_scalar_type int_poly(const base_poly &P) const
Evaluate the integral of the polynomial P on the reference element.
bgeot::pconvex_structure structure(void) const
{Structure of convex of reference.
long_scalar_type int_poly_on_face(const base_poly &P, short_type f) const
Evaluate the integral of the polynomial P on the face f of the reference element.
dim_type dim(void) const
Dimension of convex of reference.
defines and typedefs for namespace getfem
gmm::uint16_type short_type
used as the common short type integer in the library
std::shared_ptr< const convex_structure > pconvex_structure
Pointer on a convex structure description.
size_t size_type
used as the common size type in the library
pconvex_structure basic_structure(pconvex_structure cv)
Original structure (if concerned)
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.
pintegration_method exact_parallelepiped_im(size_type n)
return IM_EXACT_PARALLELEPIPED(n)
std::string name_of_int_method(pintegration_method p)
Get the string name of an integration method .
pintegration_method exact_simplex_im(size_type n)
return IM_EXACT_SIMPLEX(n)
pintegration_method exact_prism_im(size_type n)
return IM_EXACT_PRISM(n)
pintegration_method classical_approx_im(bgeot::pgeometric_trans pgt, dim_type degree)
try to find an approximate integration method for the geometric transformation pgt which is able to i...
pintegration_method classical_exact_im(bgeot::pgeometric_trans pgt)
return an exact integration method for convex type handled by pgt.
pintegration_method exact_classical_im(bgeot::pgeometric_trans pgt) IS_DEPRECATED
use classical_exact_im instead.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .
integration_method_type
the list of main integration method types
pintegration_method im_none(void)
return IM_NONE