where x component of vector_RHS is -(grad_psi)_y / x and y component is (grad_psi) / x (here psi is the another calculated scalar solution)

So, I refer to "class AdvectionField : public TensorFunction<1,dim>" in step-9.cc to form the vector_RHS

You're approaching this the wrong way. The Function<dim> and TensorFunction<dim> classes are intended to represent functions that only depend on x,y or x,y,z. But in your case, your right hand side depends on the solution itself.

I think you should look into how step-15 builds a right hand side that depends on the solution. It will probably be something like this:

  FEValuesExtractors::Scalar scalar_solution(0);
  std::vector<Tensor<1,dim>> solution_gradients (n_q_points);
  for (cell=...)
    {
       fe_values.reinit(cell);
       fe_values[scalar_solution].get_function_gradients (solution,
                                                          solution_gradients);
       for (q=...)
         {
            for (i=...)
              for (j=...)
                cell_matrix(i,j) += ...


            Tensor<1,dim> rhs;
            Point<dim> q_point = fe_values.quadrature_point(q);
            rhs[0] = solution_gradients[q][1] / q_point[0];
            rhs[1] = solution_gradients[q][0] / q_point[0];
            for (i=...)
              local_rhs(i) = fe_values[scalar_solution].value(i,q) *
                             rhs *
                             fe_values.JxW (q);
         }
     ...


Best
 W.

--
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 bange...@colostate.edu
                           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 dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.

Reply via email to