Hello,
I found out the solution and wanted to post here in case some one faces 
similar problem.

In my case, the right hand side was already assembled and there was no way 
that it could be assembled during the calculations as done usually.
Think of it as reading an assembled rhs (without application of 
Dirichlet bcs) from a file.

I found out, that if the constraints are homogeneous, I can use 
constraints.distribute_local_to_global as described below to assemble this 
rhs and apply bcs / constraints at the same time.

    std::map<unsigned int, bool> dof_touched;

    for (; cell != endc; ++cell)
        if (cell->is_locally_owned())
            for (unsigned int v = 0; v < GeometryInfo<dim>::
vertices_per_cell; ++v)
                for (unsigned int d = 0; d < dim; ++d)
                    dof_touched[cell->vertex_dof_index(v, d)] = false;


    cell = dof_handler.begin_active();

    const unsigned int   dofs_per_cell = fe.dofs_per_cell;

    std::vector<types::global_dof_index> local_dof_indices (dofs_per_cell);
    Vector<double> cell_rhs(dofs_per_cell);


for (; cell != endc; ++cell)
        if (cell->is_locally_owned())
        {
           
            cell_rhs = 0;
            cell->get_dof_indices (local_dof_indices);
            for (unsigned int i = 0; i < local_dof_indices.size(); ++i)
            {
                unsigned int dof = local_dof_indices[i];
               
                if (dof_touched[dof] == false && locally_owned_dofs.
is_element(dof))
                {
                    dof_touched[dof] = true;
                    cell_rhs(i) = assembled_rhs[dof];
                }
            }


            constraints.distribute_local_to_global (cell_rhs,
                                                    local_dof_indices,
                                                    system_rhs);
        }






On Monday, September 12, 2016 at 1:35:52 PM UTC-4, RAJAT ARORA wrote:
>
> Hi Daniel,
>
> Thanks for the reply.
>
> This indeed is one way of doing it. 
>
> I was wondering if there is any other way in which this can be done.
>
> Thanks again.
>
>
> On Monday, September 12, 2016 at 5:04:47 AM UTC-4, Daniel Arndt wrote:
>>
>> Rajat,
>>
>> The first equation is the usual equilibrium equation. It is solved to 
>>> give displacements. The stiffness matrix is then multiplied by this global 
>>> displacement vector to give the global traction vector which forms the rhs 
>>> for the second equation.
>>>
>> The easiest solution to deal with constraints to deal with them while 
>> assembling the matrix. Therefore, I would not compute the right-hand side 
>> for the second equation beforehand.
>> Instead just do it while assembling the second equation and use 
>> constraints.distribute_local_to_global().
>>
>> Best,
>> Daniel
>>
>

-- 
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