Hello Markus,

I implemented create_point_source_vector for nonprimitive elements.
Can you have a glance at my implementation? Whether it is correct way in general or not?
Thanks.

template <int dim>
void create_point_source_vector(const DoFHandler<dim> &dof_handler,
                                const Point<dim> &p,
                                Vector<double> &rhs_vector,
                                const double             factor,
                                const Tensor<1,dim> &orientation)
{
    std::pair<typename DoFHandler<dim>::active_cell_iterator, Point<dim> >
      cell_point =
GridTools::find_active_cell_around_point (StaticMappingQ1<dim>::mapping, dof_handler, p);

Quadrature<dim> q(GeometryInfo<dim>::project_to_unit_cell(cell_point.second));

    rhs_vector = 0;

FEValues<dim> fe_values(dof_handler.get_fe(), q, UpdateFlags(update_values));
    fe_values.reinit(cell_point.first);

    const unsigned int dofs_per_cell = dof_handler.get_fe().dofs_per_cell;

    std::vector<unsigned int> local_dof_indices (dofs_per_cell);
    cell_point.first->get_dof_indices (local_dof_indices);

    for(unsigned int i=0; i<dofs_per_cell; i++)
    {
        Tensor<1, dim> value;
        value[0] = fe_values.shape_value_component(i,0,0);
        value[1] = fe_values.shape_value_component(i,0,1);
        if (dim == 3)
            value[2] = fe_values.shape_value_component(i,0,2);

        rhs_vector(local_dof_indices[i]) = factor * orientation * value;

    }
}

On 26.04.2012 13:16, Markus Bürg wrote:
Hello Alexander,

no, you have a Vector<double> as destination object. Thus, you have to save a scalar, not a vector.

Best Regards,
Markus



On 26.04.2012 08:59, Alexander Grayver wrote:
Hallo Markus,

Orientation vector is defined based on experiment you want to reproduce.
For instance, modeling any finite length dipole of specified power can be done as:

/*f*(*p*) = a * *n* * \delta (*p* - p_s);/
*n* = {n_x n_y n_z} -- orientation vector;
a = I * L; I -- current, L -- length

Same can be done for magnetic dipoles (antennas).
So this is quite convenient and many real sources can be well approximated by this vector function.

However, I'm not sure now that the dot product of orientation vector and shape function is correct representation. The dot product gives us scalar, whereas we need vector with each component assigned to individual field components of FE_Nedelec. Is this still correct?

Thanks.


_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii


--
Regards,
Alexander

_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii

Reply via email to