37 #ifndef GETFEM_DERIVATIVES_H__
38 #define GETFEM_DERIVATIVES_H__
60 template<
class VECT1,
class VECT2>
62 const VECT1 &UU, VECT2 &VV) {
63 typedef typename gmm::linalg_traits<VECT1>::value_type T;
68 size_type qqdimt = qdim * N / target_qdim;
72 mf.extend_vector(UU, U);
75 "meshes are different.");
76 GMM_ASSERT1(target_qdim == qdim*N || target_qdim == qdim
77 || target_qdim == 1,
"invalid Qdim for gradient mesh_fem");
82 bgeot::pgeotrans_precomp pgp = NULL;
83 pfem_precomp pfp = NULL;
84 pfem pf, pf_target, pf_old = NULL, pf_targetold = NULL;
87 for (dal::bv_visitor cv(mf_target.
convex_index()); !cv.finished(); ++cv) {
92 GMM_ASSERT1(!(pf_target->need_G()) && pf_target->is_lagrange(),
93 "finite element target not convenient");
95 bgeot::vectors_to_base_matrix(G, mf.
linked_mesh().points_of_convex(cv));
98 if (pf_targetold != pf_target) {
99 pgp = bgeot::geotrans_precomp(pgt, pf_target->node_tab(cv), pf_target);
101 pf_targetold = pf_target;
104 pfp =
fem_precomp(pf, pf_target->node_tab(cv), pf_target);
108 gmm::dense_matrix<T> grad(N,qdim), gradt(qdim,N);
114 for (
size_type j = 0; j < pf_target->nb_dof(cv); ++j) {
118 pf->interpolation_grad(ctx, coeff, gradt, dim_type(qdim));
119 gmm::copy(gmm::transposed(gradt),grad);
120 std::copy(grad.begin(), grad.end(), V.begin() + dof_t);
124 mf_target.reduce_vector(V, VV);
142 template<
class VECT1,
class VECT2>
144 const VECT1 &UU, VECT2 &VV) {
145 typedef typename gmm::linalg_traits<VECT1>::value_type T;
150 size_type qqdimt = qdim * N * N / target_qdim;
154 mf.extend_vector(UU, U);
157 "meshes are different.");
158 GMM_ASSERT1(target_qdim == qdim || target_qdim == 1,
159 "invalid Qdim for gradient mesh_fem");
161 std::vector<T> coeff;
163 bgeot::pgeotrans_precomp pgp = NULL;
164 pfem_precomp pfp = NULL;
165 pfem pf, pf_target, pf_old = NULL, pf_targetold = NULL;
168 for (dal::bv_visitor cv(mf_target.
convex_index()); !cv.finished(); ++cv) {
171 GMM_ASSERT1(!(pf_target->need_G()) && pf_target->is_lagrange(),
172 "finite element target not convenient");
174 bgeot::vectors_to_base_matrix(G, mf.
linked_mesh().points_of_convex(cv));
177 if (pf_targetold != pf_target) {
178 pgp = bgeot::geotrans_precomp(pgt, pf_target->node_tab(cv), pf_target);
180 pf_targetold = pf_target;
183 pfp =
fem_precomp(pf, pf_target->node_tab(cv), pf_target);
187 gmm::dense_matrix<T> hess(N*N,qdim), hesst(qdim,N*N);
193 for (
size_type j = 0; j < pf_target->nb_dof(cv); ++j) {
197 pf->interpolation_hess(ctx, coeff, hesst, dim_type(qdim));
198 gmm::copy(gmm::transposed(hesst), hess);
199 std::copy(hess.begin(), hess.end(), V.begin() + dof_t);
203 mf_target.reduce_vector(V, VV);
209 template <
typename VEC1,
typename VEC2>
212 const VEC1 &U, VEC2 &VM,
214 dal::bit_vector bv; bv.add(0, mf_vm.
nb_dof());
218 template <
typename VEC1,
typename VEC2>
221 const VEC1 &U, VEC2 &VM,
222 const dal::bit_vector &mf_vm_dofs,
227 std::vector<scalar_type> DU(mf_vm.
nb_dof() * N * N);
231 GMM_ASSERT1(!mf_vm.
is_reduced(),
"Sorry, to be done");
233 scalar_type vm_min = 0., vm_max = 0.;
234 for (dal::bv_visitor i(mf_vm_dofs); !i.finished(); ++i) {
236 scalar_type sdiag = 0.;
237 for (
unsigned j=0; j < N; ++j) {
238 sdiag += DU[i*N*N + j*N + j];
239 for (
unsigned k=0; k < N; ++k) {
240 scalar_type e = .5*(DU[i*N*N + j*N + k] + DU[i*N*N + k*N + j]);
244 VM[i] -= 1./N * sdiag * sdiag;
245 vm_min = (i == 0 ? VM[0] : std::min(vm_min, VM[i]));
246 vm_max = (i == 0 ? VM[0] : std::max(vm_max, VM[i]));
248 cout <<
"Von Mises : min=" << 4*mu*mu*vm_min <<
", max="
249 << 4*mu*mu*vm_max <<
"\n";
250 gmm::scale(VM, 4*mu*mu);
257 template <
typename VEC1,
typename VEC2,
typename VEC3>
260 const VEC1 &U, VEC2 &VM,
267 typedef typename gmm::linalg_traits<VEC1>::value_type T;
269 std::vector<T> GRAD(mf_vm.
nb_dof()*N*N),
271 base_matrix sigma(N,N);
277 GMM_ASSERT1(!mf_vm.
is_reduced(),
"Sorry, to be done");
278 GMM_ASSERT1(N>=2 && N<=3,
"Only for 2D and 3D");
281 scalar_type trE = 0, diag = 0;
282 for (
unsigned j = 0; j < N; ++j)
283 trE += GRAD[i*N*N + j*N + j];
285 diag = LAMBDA[i]*trE;
287 diag = (-2./3.)*MU[i]*trE;
288 for (
unsigned j = 0; j < N; ++j) {
289 for (
unsigned k = 0; k < N; ++k) {
290 scalar_type eps = (GRAD[i*N*N + j*N + k] +
291 GRAD[i*N*N + k*N + j]);
292 sigma(j,k) = MU[i]*eps;
300 VM[i] = sqrt((3./2.)*gmm::mat_euclidean_norm_sqr(sigma));
302 VM[i] = sqrt((3./2.)*(gmm::mat_euclidean_norm_sqr(sigma) + diag*diag));
305 gmm::symmetric_qr_algorithm(sigma, eig);
306 std::sort(eig.begin(), eig.end());
307 VM[i] = eig.back() - eig.front();
void set_ii(size_type ii__)
change the current point (assuming a geotrans_precomp_ is used)
structure passed as the argument of fem interpolation functions.
Describe a finite element method linked to a mesh.
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
virtual dim_type get_qdim() const
Return the Q dimension.
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
virtual size_type nb_basic_dof() const
Return the total number of basic degrees of freedom (before the optional reduction).
virtual pfem fem_of_element(size_type cv) const
Return the basic fem associated with an element (if no fem is associated, the function will crash!...
const dal::bit_vector & convex_index() const
Get the set of convexes where a finite element has been assigned.
bool is_reduced() const
Return true if a reduction matrix is applied to the dofs.
Interpolation of fields from a mesh_fem onto another.
Define the getfem::mesh_fem class.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
pfem_precomp fem_precomp(pfem pf, bgeot::pstored_point_tab pspt, dal::pstatic_stored_object dep)
Handles precomputations for FEM.
gmm::uint16_type short_type
used as the common short type integer in the library
size_t size_type
used as the common size type in the library
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
GEneric Tool for Finite Element Methods.
void compute_gradient(const mesh_fem &mf, const mesh_fem &mf_target, const VECT1 &UU, VECT2 &VV)
Compute the gradient of a field on a getfem::mesh_fem.
void interpolation_von_mises_or_tresca(const getfem::mesh_fem &mf_u, const getfem::mesh_fem &mf_vm, const VEC1 &U, VEC2 &VM, const getfem::mesh_fem &mf_lambda, const VEC3 &lambda, const getfem::mesh_fem &mf_mu, const VEC3 &mu, bool tresca)
Compute the Von-Mises stress of a field (valid for linearized elasticity in 2D and 3D)
void interpolation_von_mises(const getfem::mesh_fem &mf_u, const getfem::mesh_fem &mf_vm, const VEC1 &U, VEC2 &VM, scalar_type mu=1)
Compute the Von-Mises stress of a field (only valid for linearized elasticity in 3D)
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
void slice_vector_on_basic_dof_of_element(const mesh_fem &mf, const VEC1 &vec, size_type cv, VEC2 &coeff, size_type qmult1=size_type(-1), size_type qmult2=size_type(-1))
Given a mesh_fem.
void compute_hessian(const mesh_fem &mf, const mesh_fem &mf_target, const VECT1 &UU, VECT2 &VV)
Compute the hessian of a field on a getfem::mesh_fem.