On 4/14/10 10:01 AM, Roy Stogner wrote:
>
> On Tue, 13 Apr 2010, Boyce Griffith wrote:
>
>> system.solution->init(system.solution->size(),
>> system.solution->local_size(),
>> ghost_dofs, GHOSTED);
>
> Why create your own ghost_dofs variable rather than just passing in
> system.get_dof_map().get_send_list()? How do you construct
> ghost_dofs?

I am not sure I understand exactly which DOFs wind up in get_send_list().

>> [0] /Users/griffith/sfw/libmesh/include/numerics/petsc_vector.h, line
>> 787, compiled Apr 13 2010 at 16:41:48
>>
>> libmesh_assert(n_local == 0 || n_local == n || !ghost.empty());
>>
>> If I comment out this assertion, everything seems to work OK.
>>
>> Am I doing something wrong here?
>
> Possibly. If you're passing in an empty ghost_dofs in a non-serial
> computation, that's almost certainly a mistake - under what
> conditions, other than the degenerate cases in the assertion, does one
> processor really need *no* data from any other?
>
> It could also be we've got an overzealous assertion, if there's a
> possibility I missed.

What I am trying to do is to compute directly the interpolation of the 
solution at the quadrature points.

I am probably being a moron, but what I'm doing is something like:

     NumericVector<Real> V;
     MeshBase::const_element_iterator el =
         mesh.active_local_elements_begin();
     const MeshBase::const_element_iterator end_el =
         mesh.active_local_elements_end();
     for ( ; el != end_el; ++el)
     {
         const Elem* const elem = *el;
         fe->reinit(elem);
         dof_map.dof_indices(elem, dof_indices);
         for (unsigned int qp = 0; qp < qrule.n_points(); ++qp)
         {
             double V_qp = 0.0;
             for (unsigned int i = 0; i < phi.size(); ++i)
             {
                 V_qp += V(dof_indices[i])*phi[i][qp];
             }
             // do stuff with V_qp
         }
     }

For this to work, obviously I need to have access to the nodal values at 
each node in each local element.

So what I do is to setup the ghost DOFs to be whatever global DOFs are 
needed on each processor which are not already local to that processor.

If all of the global DOFs associated with all of the local elements are 
local to the processor, then the vector ghost DOFs associated with that 
process will be empty.  This seems like a plausible use case to me, but 
results in this assertion failure.

Is there a better way to go about this?

Thanks,

-- Boyce

------------------------------------------------------------------------------
Download Intel&#174; Parallel Studio Eval
Try the new software tools for yourself. Speed compiling, find bugs
proactively, and fine-tune applications for parallel performance.
See why Intel Parallel Studio got high marks during beta.
http://p.sf.net/sfu/intel-sw-dev
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users

Reply via email to