Accessing raw data in MeshWorker?

78 views
Skip to first unread message

Scott Miller

unread,
Oct 25, 2012, 3:18:05 PM10/25/12
to dea...@googlegroups.com
Hi all,

I am using MeshWorker for some explicit RKDG code.  I proceed using 2 MeshWorker::integrate loops, both with Assembler::ResidualSimple<Vector<double> >.  The first integration covers all faces in my triangulation.  The second loop computes the volume contributions, including a mass matrix which I can then invert locally and update my "solution" vector.  This seems to work out very well in that I do not need any global matrices/sparsity patterns/etc.  But, I am having problems in my implementation within MeshWorker:  On each element (in my volume integration), I need access to the local component of the previously computed/assembled residual vector, e.g. something that in a non-MW code would look like "residual(local_cell_dofs[i])".  Is there currently any way to access the raw data from the IntegrationInfo or DoFInfo constructs?  I was unable to find any interface for doing so.  If nothing exists, should we add a raw data accessor function that take either an index or name as the argument?

Thanks.

-Scott

Bärbel Janssen

unread,
Oct 25, 2012, 3:30:38 PM10/25/12
to dea...@googlegroups.com
Hi Scott,

you can add your old residual to the NamedData that you pass to your
loop and then use ist within the MeshWorker. Does that help?

Best,
Bärbel

2012/10/25 Scott Miller <miller....@gmail.com>:
> --
> 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
>
>

Scott Miller

unread,
Oct 25, 2012, 3:36:47 PM10/25/12
to dea...@googlegroups.com
Hi Barbel,

I have added it to the NamedData that I give to MeshWorker::IntegrationInfoBox. But how do I access that from within when I'm on a cell? I can get the solution values, gradients, etc at the quadrature points fine, but I'm not sure how to get the raw residual vector.

-Scott

Bärbel Janssen

unread,
Oct 25, 2012, 4:23:23 PM10/25/12
to dea...@googlegroups.com
Hi Scott,

let's assume that your NamedData is called data. Then you can do
data(data.find("solution")), or you can print it like
data.print(std::cout) to see what was added to the NamedData. This way
you can access all the data you added. Does that work?

Scott Miller

unread,
Oct 25, 2012, 5:01:02 PM10/25/12
to dea...@googlegroups.com
Hi Barbel,

Let us take a step back, and I ask: How can I access the NamedData object from a MeshWorker::IntegrationInfo object?

Thanks.

-Scott

Bärbel Janssen

unread,
Oct 25, 2012, 5:19:32 PM10/25/12
to dea...@googlegroups.com
Hi Scott,

you can get it by info.values[position of data in NamedData] and the
same with gradients and second derivatives.

Scott T. Miller

unread,
Oct 25, 2012, 5:47:23 PM10/25/12
to dea...@googlegroups.com, Bärbel Janssen
Hi,

I thought the info.values accessed the values of the finite element solution evaluated at the quadrature points, not the "raw" data directly.

-Scott

Guido Kanschat

unread,
Oct 26, 2012, 6:07:16 AM10/26/12
to dea...@googlegroups.com, Bärbel Janssen
Dear Scott,

I have written a class VectorInfo analog to IntegrationInfo which does
what you want. I attached the file and am also adding it to the
subversion archive.

Best,
Guido
vector_info.h

Scott T. Miller

unread,
Oct 26, 2012, 11:25:23 AM10/26/12
to dea...@googlegroups.com
Hi Guido,

Many thanks! I will test it out.

-Scott
> --
> 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
>
>
> <vector_info.h>

Scott Miller

unread,
Oct 30, 2012, 9:56:10 AM10/30/12
to dea...@googlegroups.com, Bärbel Janssen
Dear Guido,

Sorry for the delay in my testing your new class. I am having some difficulty in understanding its usage. Where would I make and/or use a VectorInfoBox? I still need the IntegrationInfo and DoFInfo in my assemble function. I also see that you have an access function ::operator(unsigned int i). If I have a VectorInfo object vec_info, and I want to do "Vector<double> my_vec = vec_info(i).block(0)", what would I be passing as the value i?

Thanks again,

-Scott

On Oct 26, 2012, at 6:07 AM, Guido Kanschat wrote:

> --
> 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
>
>
> <vector_info.h>

Guido Kanschat

unread,
Oct 30, 2012, 10:03:14 AM10/30/12
to dea...@googlegroups.com, Bärbel Janssen

Oops, maybe I should have added a little information on its use.

It is a replacement for IntegrationInfo. Or, VectorInfoBox for IntegrationInfoBox.

Do you need the functionality of both?

Best,
Guido

Scott Miller

unread,
Oct 30, 2012, 10:06:37 AM10/30/12
to dea...@googlegroups.com, Bärbel Janssen
Hi Guido,

Yes, I do need the functionality of both.  I have a non-linear Jacobian that I linearize and invert on the interior of each element, and the multiply the cell vector coming from the residual of the face terms.

Thanks.

-Scott
Reply all
Reply to author
Forward
0 new messages