On 7/15/23 11:33, Mohammad Amir Kiani Fordoei wrote:

const Tensor<1,dim> n = fe_face_values.normal_vector(q_point);
cell_rhs (i)+=( (n*traction_component*fe_face_values.shape_value(i, q_point))*
fe_face_values.JxW(q_point));
}//for q_point
}//for dof_per_cell

But in computation of cell rhs i get this error:

/home/amir/eclipse-workspace/mech1/mech.cc:983:22: error: no match for ‘operator+=’ (operand types are ‘double’ and ‘dealii::Tensor<1, 3>’)   983 |          cell_rhs (i)+=( (n*traction_component*fe_face_values.shape_value(i, q_point))*       |  ~~~~~~~~~~~~^~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
   984 |                   fe_face_values.JxW(q_point));

As error says, the normal vector with 3 component cannot be multiplied by double values in other terms.

No, that's not actually what the error message says :-) It doesn't say anything about multiplication. If you read it carefully, it says that you cannot assign a right hand side of type Tensor<1,3> to a left hand side of type double. That makes sense.

The types you have are because your right hand side is the product of
  n  -- a tensor
  traction_component -- not sure, but probably a scalar
  fe_face_values.shape_value(i, q_point))  -- a scalar, see the function's
                                              documentation
  fe_face_values.JxW(q_point) -- a scalar

What you *want* is a vector-valued shape function. You get that by "extracting" this from the fe_values object in the same way as is discussed in the vector-valued documentation module. I think something like this should work:

FEValuesExtractors::Vector displacement(0);
cell_rhs (i)+=( (n*traction_component*fe_face_values[displacement].shape_value(i, q_point))*
fe_face_values.JxW(q_point));

Best
 W.

--
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                           www: http://www.math.colostate.edu/~bangerth/


--
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 [email protected].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/9f48ee6b-cf61-499b-985e-268524743ee0%40colostate.edu.

Reply via email to