> On 17 Feb 2017, at 14:46, Franco Milicchio <franco.milicc...@gmail.com> wrote:
> 
> I agree I have a translation problem here. In practice, I am integrating a 
> weak problem as this (in Fenics' UFL terms, but easy to recognize):
> 
> 
>     inner( sigma(u), eps(v) ) * dx = inner(f, v) * dx + inner(f1, v) * ds(1)


I’m guessing that 

f1 = sigma(u) n 

is the Neumann data on the boundary. Then if you have f1 as a dim-dimensional 
vector field, it already contains n. 

the translation from ufl to deal.II is easy:


ds(1) =>  fe_face_values.JxW(f_q_point); 
v => fe_face_values.shape_value(i, f_q_point)
f1 =>  rhs_face_values[f_q_point]

since fe_face_values.shape_value(i, f_q_point) refers to the component_i of v, 
then

inner(f1, v) *ds(1) = sum over q of 

fe_face_values.shape_value(i, f_q_point) * 
rhs_face_values[f_q_point][component_i]) * fe_face_values.JxW(f_q_point); 


> where f is a bulk force (zero in my code), f1 is a (distributed) force 
> applied to the top boundary. The dx and ds terms refer to the integral on the 
> domain and boundary; u is the displacement, and v is the test function. Sigma 
> and eps are functions that compute, obviously, stress and strain tensors.
> 
> As I said, Fenics hid everything as you see above, but I am very eager to let 
> it go as soon as possible.
> 
>  
> and while you can certainly come up with good reasons for this to work out, 
> I’m unsure about what you want to achieve. If your boundary conditions are 
> 
> ((grad(u)).n,v) = u_i,j n_j v_i 
> 
> and you impose 
> 
> (grad(u).n)_i = f_i 
> 
> on the boundary, then you should have 
> 
> cell_rhs(i) += (    fe_face_values.shape_value(i, f_q_point) * 
>                     rhs_face_values[f_q_point][component_i]) * 
> fe_face_values.JxW(f_q_point); 
> 
> 
> My question here is this: why am I not using the normal to the face?

You are not using it, because you are imposing

sigma(u)n = f

on the boundary (where f is a dim-dimensional vector field on the boundary). 
Even in your FENICS implementation you are not using the normal… 

inner(f1,v)*ds(1) 

does not contain the normal anywhere, only the integral of f_i v_i on the faces 
on the boundary (I guess with id 1).


By the way this is the same for the Poisson problem. If you impose 

u,i n_i = f

then your integral on the boundary becomes

f*v*ds(1)












-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to