Dear Prof. Bangerth,
                                 Thank you for the guidance and correction. 
I got the point now that for Neuman conditions I shall iterate for face 
cells along desired edge for imposing condition like in the previous case 
of step-8. I also went through the step-7, thank you for the guidance it 
made me understand more about the face_cell_values for making rhs. As far 
as I have understood the step-7 is iterating through each cell and 
implementing simply the Neuman condition at the face cells:

for (; cell!=endc; ++cell)
{
cell_matrix = 0;
cell_rhs = 0;
fe_values.reinit 
<https://www.dealii.org/current/doxygen/deal.II/classFEValues.html#aec8f5b8b3e4c5dcf35dfd029a1ecbbd0>
 
(cell);
right_hand_side.value_list (fe_values.get_quadrature_points 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#a5f8732ebe2d3c6746f6de26a79cb1e45>
(),
rhs_values);
for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
for (unsigned int i=0; i<dofs_per_cell; ++i)
{
for (unsigned int j=0; j<dofs_per_cell; ++j)
cell_matrix(i,j) += ((fe_values.shape_grad 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ac4cee7628c2903a89c5c399fddeb00a5>(i,q_point)
 
*
fe_values.shape_grad 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ac4cee7628c2903a89c5c399fddeb00a5>
(j,q_point)
+
fe_values.shape_value 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#abe4de48ff59778bb82a0ec13037804aa>(i,q_point)
 
*
fe_values.shape_value 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#abe4de48ff59778bb82a0ec13037804aa>(j,q_point))
 
*
fe_values.JxW 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ad097580a2f71878695096cc73b271b9d>
(q_point));
cell_rhs(i) += (fe_values.shape_value 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#abe4de48ff59778bb82a0ec13037804aa>(i,q_point)
 
*
rhs_values [q_point] *
fe_values.JxW 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ad097580a2f71878695096cc73b271b9d>
(q_point));
}
for (unsigned int face_number=0; 
face_number<GeometryInfo<dim>::faces_per_cell; ++face_number)
if (cell->face(face_number)->at_boundary()
&&
(cell->face(face_number)->boundary_id() == 1))
{
fe_face_values.reinit 
<https://www.dealii.org/current/doxygen/deal.II/classFEFaceValues.html#af6e079ca7429d54433343d50bd334c3c>
 
(cell, face_number);
for (unsigned int q_point=0; q_point<n_face_q_points; ++q_point)
{
const double neumann_value
= (exact_solution.gradient (fe_face_values.quadrature_point 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ab1aa3aa2940125b47df95ff82ad733b9>(q_point))
 
*
fe_face_values.normal_vector 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#a130eea0fa89263d93b20521addc830c7>
(q_point));
for (unsigned int i=0; i<dofs_per_cell; ++i)
cell_rhs(i) += (neumann_value *
fe_face_values.shape_value 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#abe4de48ff59778bb82a0ec13037804aa>(i,q_point)
 
*
fe_face_values.JxW 
<https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#ad097580a2f71878695096cc73b271b9d>
(q_point));
}
}
cell->get_dof_indices (local_dof_indices);
for (unsigned int i=0; i<dofs_per_cell; ++i)
{
for (unsigned int j=0; j<dofs_per_cell; ++j)
system_matrix.add (local_dof_indices[i],
local_dof_indices[j],
cell_matrix(i,j));
system_rhs(local_dof_indices[i]) += cell_rhs(i);
}
}


*But the step-26 which I am now concerned about is directly forming rhs 
somehow using the function:*

        RightHandSide<dim> rhs_function;
        rhs_function.set_time(time);
        VectorTools::create_right_hand_side(dof_handler,
                                            QGauss<dim>(fe.degree+1),
                                            rhs_function,
                                            tmp);
So how can I possibly implement heat flux on cell faces in this example 
case step-26 ? If I can, then I guess my problem is solved to learn 
implementing Neuman BC for Heat equation in deal.ii. Thank you very much in 
advance? 



On Wednesday, April 10, 2019 at 5:56:46 AM UTC+2, Wolfgang Bangerth wrote:
>
> On 4/9/19 10:01 AM, Muhammad Mashhood wrote: 
> > I want to apply the heat flux BC on one of the edge instead i.e. *dU(t, 
> x, 
> > y=1)/dy* at top. For this so far I could only find a method to make a 
> heat 
> > source using the already present function as follows: 
>
> I don't know how you use the class you implement, but I would try to avoid 
> all 
> confusion by naming the class RightHandSide: This is the name we generally 
> use 
> for forcing terms of the differential equation, and not for the right hand 
> side of boundary conditions. 
>
> Just like for your previous question, it is not correct to "convert" the 
> right 
> hand side of a boundary condition into a right hand side of the equation. 
> This 
> is mathematically not correct and will not lead to the correct solution. 
>
> Have you looked at step-7? This program deals with Neumann boundary 
> conditions 
> and may have exactly what you are looking for. 
>
> Best 
>   W. 
>
> -- 
> ------------------------------------------------------------------------ 
> Wolfgang Bangerth          email:                 [email protected] 
> <javascript:> 
>                             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].
For more options, visit https://groups.google.com/d/optout.

Reply via email to