Hi Giovani,
If you want to retrieve the coordinates of the quadrature points of a
convex cv, you can do:
/* get the approximate integration method for the convex */
getfem::papprox_integration pai
= mim.int_method_of_element(cv)->approx_method();
/* prepare the geotrans_interpolation_context to apply the geometric
transformation */
bgeot::base_matrix G;
bgeot::vectors_to_base_matrix(G, m.points_of_convex(cv));
bgeot::geotrans_interpolation_context c(m.trans_of_convex(cv),
pai->point(0), G);
/* apply the geotrans to each point of the integration method */
for (size_type j = 0; j < pai->nb_points_on_convex(); ++j) {
c.set_xref(pai->point(j));
std::cout << " integration node " << j << " for convex " << cv <<
" is " << c.xreal() << "\n";
}
Best regards,
Julien
Giovani wrote:
Good day gentlemen!
Suppose I have a getfem::mesh object that has a getfem::mesh_fem
linked to it, and also a getfem::mesh_im that's also linked.
I would like to build a new getfem::mesh, whose nodes are located at
the quadrature points of the previous mesh. Is there any simple way of
retrieving a list with the global coordinates of those quadrature points?
Best Regards,
Giovani M. Faccin
_______________________________________________
Getfem-users mailing list
[email protected]
https://mail.gna.org/listinfo/getfem-users
--
Julien Pommier, bureau 111
GMM, INSA Toulouse, tél:05 61 55 93 42
_______________________________________________
Getfem-users mailing list
[email protected]
https://mail.gna.org/listinfo/getfem-users