[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]
Re: [Getfem-users] Build a mesh with quadrature points
From: |
julien pommier |
Subject: |
Re: [Getfem-users] Build a mesh with quadrature points |
Date: |
Mon, 19 Mar 2007 11:38:58 +0100 |
User-agent: |
Thunderbird 1.5.0.10 (X11/20070307) |
Hi Giovani,
You can interpolate a mesh_fem at arbitrary locations using the function
void getfem::interpolation(const mesh_fem &mf_source, mesh_trans_inv &mti,
const VECTU &U, VECTV &V, bool extrapolation = false)
It is defined in getfem_interpolation.h
Its mesh_trans_inv argument is where you will store the list of global
coordinates to interpolate. U is the data on mf_source (U.size() ==
mf_source.nb_dof()), V is the interpolated data on the points stored in
mti (V.size() == mti.nb_points()*mf_source.get_qdim()).
mti is filled like that:
mesh_trans_inv mti(mf_source.linked_mesh());
mti.add_point(base_node(42,42));
etc..
Now if you want to be able to perform this interpolation very quickly
for a set of data vectors U, you can use the function
getfem::interpolation(const mesh_fem &mf_source,
mesh_trans_inv &mti,
const VECTU &U, VECTV &V, MAT &M,
int version, bool extrapolation = false)
which is very similar, except that it, if you set version to 1, instead
of filling V, it will fill the (preferably sparse) matrix M in such way
that each row of M will contain the value of each base function on a
point of mti. So this is indeed the answer to your first question. And
of course, when you have that matrix M, you just have to multiply it
with a vector U to obtain the interpolated values V, so you can perform
many interpolations very quickly.
I hope I was clear ;)
Best regards,
Julien
Giovani wrote:
Julien
Thank you very much! That piece of code is helping a lot!
If I may ask, I have another issue that requires enlightenment:
Suppose you have a getfem::mesh, a getfem::mesh_fem linked with the
mesh, and a point in global coordinates whose location is not
identical to any of the mesh nodes. This point is, however, inside the
object that the mesh represents, so I can interpolate data from the
mesh nodes to it.
--> How can I evaluate the shape functions of a single, specific
getfem::mesh_fem dof at this point's coordinates?
--> If I multiply the result from the previous question by the value
of a field prescribed at that dof, I'll then have the contribution of
that dof to the interpolation of the field at the point. Just to be
sure, this result will be normalized, right? I mean, if I repeat this
for all dofs and simply sum the results, I'd obtain the same
interpolation that I'd get if I used getfem::interpolation instead,
although this explicit sum is probably more computationaly expensive
than using getfem::interpolation. Is that it?
Best regards,
Giovani
julien pommier escreveu:
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
address@hidden
https://mail.gna.org/listinfo/getfem-users
_______________________________________________
Getfem-users mailing list
address@hidden
https://mail.gna.org/listinfo/getfem-users
--
Julien Pommier, bureau 111
GMM, INSA Toulouse, tél:05 61 55 93 42