Looping over cells and faces: MeshWorker vs WorkStream

71 views
Skip to first unread message

Praveen C

unread,
May 31, 2021, 12:21:33 AM5/31/21
to Deal. II Googlegroup
Dear all

I have a DG-type scheme which is totally quadrature-free and matrix-free. I need to first loop over cells and then loop over faces. I know MeshWorker can be used for this, but it seems to be not used by many people and I suspect it may not be maintained well.

Is it possible to use WorkStream to loop over faces ?

Thanks
praveen

Jean-Paul Pelteret

unread,
May 31, 2021, 1:57:49 PM5/31/21
to dea...@googlegroups.com
Hi Praveen,

Is it possible to use WorkStream to loop over faces ?

It is, in the sense that you can use it to multithread the loop over all cells, and then loop over all all faces on a cell yourself. 

This has been automated in the mesh_loop() function (in fact, mesh_loop() is used as a part of the backend to MeshWorker). There you can specify a boundary worker that will be called either before all after the cells — you get to choose using the MeshWorker::AssembleFlags . It’s quite a convenient assembly device, and its being put to use in most of the new tutorials. There are a few tutorials that use mesh_loop() already. Just look in the keywords of the tutorial list to find them. I’m sure that there must also be a DG example in the test suite somewhere, if you need it.

I know MeshWorker can be used for this, but it seems to be not used by many people and I suspect it may not be maintained well. 

Yes, it seems not to be used terribly often (at least, the questions on the mailing list are infrequent, and the number of people that can answer those questions seems somewhat limited). I advocate using mesh_loop() in preference to MeshWorker.

Best,
Jean-Paul

--
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+un...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/98B64908-C102-4C4F-B18E-99066E32C1A4%40gmail.com.

Jean-Paul Pelteret

unread,
May 31, 2021, 2:54:08 PM5/31/21
to dea...@googlegroups.com
Hi Praveen,

Correction, I meant to say face worker instead of boundary worker. Step-74 is a DG code that uses mesh_loop() with a face worker, and also makes use of the nice FEInterfaceValues class that has also very recently been extended to use FEValuesExtractors (see item 17 in the specific improvements section for the upcoming release). 

Best,
Jean-Paul

Praveen C

unread,
May 31, 2021, 11:36:39 PM5/31/21
to Deal. II Googlegroup
Thanks a lot for these hints.

I have used MeshWorker::loop many times but not the newer MeshWorker::mesh_loop.

In mesh_loop, on an interior face, I need to compute rhs vector for the two adjacent cells. 

I dont have any matrix since I am doing explicit DG scheme.

To assemble rhs vector on interior face, I am guessing I have to store it like this

struct CopyDataFace
{
};
struct CopyData
{
Vector<double> cell_rhs;
std::vector<CopyDataFace> face_data;

};

It this how it should be used ?

Thanks
praveen

Wolfgang Bangerth

unread,
Jun 1, 2021, 12:06:42 AM6/1/21
to dea...@googlegroups.com, Praveen C
On 5/31/21 9:36 PM, Praveen C wrote:
>
> It this how it should be used ?
>

Something like this. You might want to look at step-47 which assembles face
terms for interior faces.

Best
W.

--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Corbin Foucart

unread,
Jul 19, 2021, 3:50:27 PM7/19/21
to deal.II User Group
I'm writing a similar code to Praveen, and had some additional questions, having looked at the documentation for mesh_loop as well as step-47, step-51, and the code gallery entry for the second-order LDG solver.

For a scheme that has a mixed-formulation element-local system (Q,U) as does the LDG scheme in the code gallery, it is desirable to instead eliminate the gradient Q on each cell and assemble a linear system only in U. This amounts to a Schur complement taken on each cell--that is,
[[ A       B ]]  [Q]  = [ 0 ]
[[ -B^T  C ]]  [U]     [ F ]

and a cell contribution to the LHS of the global linear system KU = L, where K = C + B^T A^{-1} B. My question is the following:
  • It seems that since cell interior data (A inverse) is needed to premultiply the B operator (which contains face terms), it makes more sense to loop over each cell and handle the faces of the cell manually, like in step-51, which performs a similar element-local static condensation.
  • In the spirit of step-51, my idea is to use mesh_loop, but only use the cell-assembly part and pass dummy functions to the face_worker and boundary_worker arguments. Is that a sensible way to do it, or should I be using a different MeshWorker function?
  • Or alternatively, is there a way to use the face_worker function, but preserve the data on the cell to the left and right? It seems expensive to recompute the inverse mass matrix on every interface, for each adjacent cell.
Best,
Corbin

Wolfgang Bangerth

unread,
Jul 19, 2021, 5:24:51 PM7/19/21
to dea...@googlegroups.com
On 7/19/21 1:50 PM, Corbin Foucart wrote:
> * It seems that since cell interior data (A inverse) is needed to
> premultiply the B operator (which contains face terms), it makes more
> sense to loop over each cell and handle the faces of the cell manually,
> like in step-51, which performs a similar element-local static condensation.
> * In the spirit of step-51, my idea is to use mesh_loop, but only use the
> cell-assembly part and pass dummy functions to the face_worker and
> boundary_worker arguments. Is that a sensible way to do it, or should I be
> using a different MeshWorker function?
> * Or alternatively, is there a way to use the face_worker function, but
> preserve the data on the cell to the left and right? It seems expensive to
> recompute the inverse mass matrix on every interface, for each adjacent cell.

All of these seem reasonable approaches. Unless you have concrete evidence
that whatever performance you get, going with the easiest-to-implement
approach is always the right choice!

Best
Wolfgang
Reply all
Reply to author
Forward
0 new messages