<html>
  <head>
    <meta content="text/html; charset=windows-1252"
      http-equiv="Content-Type">
  </head>
  <body bgcolor="#FFFFFF" text="#000000">
    Hi Arya,<br>
    <br>
    there is a mail by Jö Fahlke from May 16 2011 from this mailing
    list, where it is explicitly explained how to obtain an index for an
    arbitrary intersection, in response to a question similar to yours.
    Let me just quote this mail:<br>
    <br>
    <div class="moz-text-plain" wrap="true" graphical-quote="true"
      style="font-family: -moz-fixed; font-size: 12px;" lang="x-unicode">
      <pre wrap="">Am Mon, 16. May 2011, 11:27:09 +0200 schrieb Benjamin Faigle:
</pre>
      <blockquote type="cite" style="color: #000000;">
        <pre wrap="">if i understood it right, the method
        grid.localIdSet().id(const EntityType &e)
should return a unique ID for every sort of entities, be it
elements, vertices or intersections. It works out pretty well with
elements, as shown for example in the Dune Grid Howto
        IdType idf = idset.id (* ep );

However, if I try the same with an *intersection, I get the error message:
...dune-grid/dune/grid/common/indexidset.hh:394:43: error:
‘codimension’ is not a member of ‘Dune::IntersectionIterator<const
Dune::UGGrid<2>, Dune::UGGridLeafIntersectionIterator,
Dune::UGGridLeafIntersection>’

I am using dune 2.0.
I am sorry if it is just an incredibly silly question, and many
thanks for your help,
</pre>
      </blockquote>
      <pre wrap="">Only entities have ids.  Intersections are not entities!

All grids have entities of codimension dim (vertices) and codimension 0
(elements), and all have intersections.  Some grids have entities of other
codimensions as well.  Even if a grid has entities of codimension 1, these are
not the same as entities.  Most obvious differences:

 * Intersections know the neighbouring cells, codim 1 entities don't.

 * Intersections have an orientation, codim 1 entities don't.

 * Codim 1 entities have indices and ids and can be used to obtain array
   indices from a mapper (even if the grid does not support codim 1 entity
   objects!)  For intersections, indices, id's and maps are not provided by
   core Dune, though it is possibly to build such things yourself.

 * In non-conforming grids you may have the following situation:

     +-------+---+---+
     |       | b | d |
     |   a   +---+---+
     |       | c | e |
     +-------+---+---+

   Elements b, c, d and e all have four codim 1 subentities and four
   corresponding intersections.  Element a also has four codim 1 subentities,
   but it has five intersections.  The codim 1 entity to the right of b is the
   same as the codim 1 entity to the left of d.  You have three codim 1
   entities to the right of a however: on covering the whole of a's right
   border, and two smaller one each covering the upper and the lower half of
   a's border, respectively.

If you know you have a conforming grid, you have a correspondance between
intersections and codim 1 entities (ignoring the orientation of the
intersections).  In that case you can usually get away by using the id of the
corresponding codim 1 entity as the id of the intersection (even if the grid
does not support codim 1 entities as independent objects.  Use

  grid.globalIdSet().subId(e, i.indexInInside(), 1);

where i is the intersection and e is the codim 0 entity i was obtained from.
If you don't have e readily available, you can use

  grid.globalIdSet().subId(*i.inside(), i.indexInInside(), 1);

instead.

Bye,
Jö.


</pre>
    </div>
    I hope this helps.<br>
    Best,<br>
    Jurgis<br>
    <br>
    <br>
    <div class="moz-cite-prefix">On 06.06.2013 14:01, Arya Fallahi
      wrote:<br>
    </div>
    <blockquote
cite="mid:CAHr884wnuLikPMHh=xHAP+LLKrohKCwvrks3-vf3iRPdyDsAWQ@mail.gmail.com"
      type="cite">
      <div dir="ltr">
        <div>
          <div>
            <div>
              <div>Dear Jurgis and Steffen,<br>
                <br>
              </div>
              Following your discussion, I thought of asking you my
              question since the problem Jurgis has mentioned is somehow
              close to mine. I have written yesterday to dune but have
              received no feedback yet. <br>
              <br>
              I am indeed not completely aware about the kind of data
              Jurgis is going to attach to intersections but what I need
              to do is as follows:<br>
              <br>
            </div>
            I am using .msh format produced by gmsh to prepare my grid
            and the tags attached to the entities. In a specific
            analysis, I need to find intersections with specific tags.
            These tags for any kind of elements are supported by .msh
            and also by dgf format and ALUGrid. What I need to do is to
            call the tags of intersections from dune. The same thing as
            what intersect->boundaryId() does is needed. However,
            boundaryId() gives the tags if and only if the intersection
            is of boundary type i.e. intersect->boundaryId() = true.
            I have read the dune tutorial in detail and also searched
            through out the class documentation but could not find how
            the tags can be retrieved if intersect->boundaryId() =
            false. Do you know if this is supported or implemented at
            all in dune? <br>
            <br>
          </div>
          Thank you in advance for your help,<br>
          <br>
        </div>
        Arya<br>
      </div>
      <div class="gmail_extra"><br>
        <br>
        <div class="gmail_quote">On Thu, Jun 6, 2013 at 1:24 PM, Jurgis
          Pods <span dir="ltr"><<a moz-do-not-send="true"
              href="mailto:jurgis.pods@iwr.uni-heidelberg.de"
              target="_blank">jurgis.pods@iwr.uni-heidelberg.de</a>></span>
          wrote:<br>
          <blockquote class="gmail_quote" style="margin:0 0 0
            .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi
            Steffen,<br>
            <br>
            thanks for the hint! I actually got it to work, it turned
            out that my problem actually was that I used
            dune-multidomaingrid the wrong way!<br>
            <br>
            [Begin Off-topic]<br>
            When defining the multidomain grid like this:<br>
            <br>
            typedef Dune::MultiDomainGrid<HostGrid,Dune::mdgrid::FewSubDomainsTraits<HostGrid::dimension,2,Dune::mdgrid::CellAndVertexCodims>
            > Grid;<br>
            <br>
            <br>
            and then calling for an IntersectionIterator 'iit' the
            following:<br>
            <br>
            gv.indexSet().subIndex(*iit->inside(),
            iit->indexInInside(), 1)<br>
            <br>
            <br>
            on the multidomain gridview, everything goes well. But when
            calling this on a *subdomain* gridview:<br>
            <br>
            subgv.indexSet().subIndex(*iit->inside(),
            iit->indexInInside(), 1)<br>
            <br>
            <br>
            the program would exit and complain with<br>
            <br>
            Dune reported error: GridError
            [invoke:/home/jpods/workspace/dev/dune-multidomaingrid/dune/grid/multidomaingrid/indexsets.hh:92]:
            codim not supported<br>
            <br>
            <br>
            Apparently, the lesson learned is to only use the subIndex()
            method on the multidomain gridview.<br>
            [End Off-topic]<br>
            <br>
            Thanks again!<span class="HOEnZb"><font color="#888888"><br>
                Jurgis</font></span>
            <div class="HOEnZb">
              <div class="h5"><br>
                <br>
                On 06.06.2013 10:13, Steffen Müthing wrote:<br>
                <blockquote class="gmail_quote" style="margin:0 0 0
                  .8ex;border-left:1px #ccc solid;padding-left:1ex">
                  Hi Jurgis,<br>
                  <br>
                  you don't need actual face entity objects to do this,
                  you just need their indices (and those should be
                  provided by<br>
                  all DUNE grids, if I understand the grid paper
                  correctly - you are fine at least, because YaspGrid
                  definitely provides<br>
                  indices for faces). If your grid does not support face
                  entities, you can still get their indices with the
                  help of the cell<br>
                  entity if you know the embedding of the face in the
                  cell (usually from Intersection::indexInInside() ).
                  You can even get<br>
                  the indices of the face's subentities. Take a look at
                  IndexSet::subIndex() for the details.<br>
                  <br>
                  That said (and slightly off-topic for this list), I
                  don't have the slightest idea if IntersectionIndexSet
                  still works. That thing<br>
                  was only ever really used by Sven's mimetic FD code,
                  and that hasn't been maintained for quite a while now.
                  YMMV<br>
                  (it's definitely broken in PDELab master, though)...<br>
                  <br>
                  Steffen<br>
                  <br>
                  Am 05.06.2013 um 14:58 schrieb Jurgis Pods:<br>
                  <br>
                  <blockquote class="gmail_quote" style="margin:0 0 0
                    .8ex;border-left:1px #ccc solid;padding-left:1ex">
                    Hi Steffen,<br>
                    <br>
                    this approach (using PDELab's IntersectionIndexSet)
                    will not work on grids which do not support codim-1
                    entities, will it? Does this mean I have to use an
                    IdSet in this case?<br>
                    <br>
                    Jurgis<br>
                    <br>
                    <br>
                    On <a moz-do-not-send="true"
                      href="tel:07.11.2012%2022" value="+41711201222"
                      target="_blank">07.11.2012 22</a>:21, Steffen
                    Müthing wrote:<br>
                    <blockquote class="gmail_quote" style="margin:0 0 0
                      .8ex;border-left:1px #ccc solid;padding-left:1ex">
                      Hi Marco,<br>
                      <br>
                      as long as your grid has a conforming macro grid
                      (which is true for every Dune grid I know),<br>
                      you can always associate the data with the face
                      belonging to the more-refined<br>
                      entity (the one with the higher level), or with
                      the face belonging to the cell with<br>
                      the lower index (if the levels are the same). To
                      get the index of that face, you can<br>
                      use the subIndex() method of the IndexSet with the
                      cell and the value of indexInInside()<br>
                      or indexInOutside() from the intersection.<br>
                      <br>
                      In PDELab, there is an IntersectionIndexSet that
                      does something like that. Perhaps you<br>
                      can simply use that, or you can use it as
                      inspiration for your own implementation… ;-)<br>
                      <br>
                      Steffen<br>
                      <br>
                      <br>
                      Am 07.11.2012 um 18:03 schrieb Marco Cisternino:<br>
                      <br>
                      <br>
                      <blockquote class="gmail_quote" style="margin:0 0
                        0 .8ex;border-left:1px #ccc
                        solid;padding-left:1ex">
                        Hi Duners,<br>
                        I need to attach data to intersections. I would
                        like to avoid redundance of data (using a mapper
                        to elements and attaching a stl container to
                        store values on its intersections) and time
                        consumption to access data (using a stl map).<br>
                        I wonder if something like a mapper is available
                        for intersections. Choosing the right layout for
                        a MCMG mapper would not work, I guess, on non
                        conforming grid, doesn't it?<br>
                        Thank you very much, again.<br>
                        <br>
                        Marco<br>
                        <br>
                        <br>
                        -- <br>
                        Marco Cisternino<br>
                        Optimad Engineering s.r.l.<br>
                        <br>
                        <a moz-do-not-send="true"
                          href="http://www.optimad.it" target="_blank">www.optimad.it</a><br>
                        <a moz-do-not-send="true"
                          href="mailto:marco.cisternino@optimad.it"
                          target="_blank">marco.cisternino@optimad.it</a><br>
                        <br>
                        <a moz-do-not-send="true"
                          href="tel:%2B3901119719782"
                          value="+3901119719782" target="_blank">+3901119719782</a><br>
                        <br>
                        <br>
                        _______________________________________________<br>
                        Dune mailing list<br>
                        <br>
                        <a moz-do-not-send="true"
                          href="mailto:Dune@dune-project.org"
                          target="_blank">Dune@dune-project.org</a><br>
                        <a moz-do-not-send="true"
                          href="http://lists.dune-project.org/mailman/listinfo/dune"
                          target="_blank">http://lists.dune-project.org/mailman/listinfo/dune</a><br>
                      </blockquote>
                      Steffen Müthing<br>
                      Universität Stuttgart<br>
                      Institut für Parallele und Verteilte Systeme<br>
                      Universitätsstr. 38<br>
                      70569 Stuttgart<br>
                      Tel: <a moz-do-not-send="true"
                        href="tel:%2B49%20711%20685%2088429"
                        value="+4971168588429" target="_blank">+49 711
                        685 88429</a><br>
                      Fax: <a moz-do-not-send="true"
                        href="tel:%2B49%20711%20685%2088340"
                        value="+4971168588340" target="_blank">+49 711
                        685 88340</a><br>
                      Email:<br>
                      <a moz-do-not-send="true"
                        href="mailto:steffen.muething@ipvs.uni-stuttgart.de"
                        target="_blank">steffen.muething@ipvs.uni-stuttgart.de</a><br>
                      <br>
                      <br>
                      <br>
                      <br>
                      <br>
                      _______________________________________________<br>
                      Dune mailing list<br>
                      <br>
                      <a moz-do-not-send="true"
                        href="mailto:Dune@dune-project.org"
                        target="_blank">Dune@dune-project.org</a><br>
                      <a moz-do-not-send="true"
                        href="http://lists.dune-project.org/mailman/listinfo/dune"
                        target="_blank">http://lists.dune-project.org/mailman/listinfo/dune</a><br>
                    </blockquote>
                  </blockquote>
                  Steffen Müthing<br>
                  Universität Stuttgart<br>
                  Institut für Parallele und Verteilte Systeme<br>
                  Universitätsstr. 38<br>
                  70569 Stuttgart<br>
                  Tel: <a moz-do-not-send="true"
                    href="tel:%2B49%20711%20685%2088429"
                    value="+4971168588429" target="_blank">+49 711 685
                    88429</a><br>
                  Fax: <a moz-do-not-send="true"
                    href="tel:%2B49%20711%20685%2088340"
                    value="+4971168588340" target="_blank">+49 711 685
                    88340</a><br>
                  Email: <a moz-do-not-send="true"
                    href="mailto:steffen.muething@ipvs.uni-stuttgart.de"
                    target="_blank">steffen.muething@ipvs.uni-stuttgart.de</a><br>
                  <br>
                </blockquote>
                <br>
                <br>
                _______________________________________________<br>
                Dune mailing list<br>
                <a moz-do-not-send="true"
                  href="mailto:Dune@dune-project.org" target="_blank">Dune@dune-project.org</a><br>
                <a moz-do-not-send="true"
                  href="http://lists.dune-project.org/mailman/listinfo/dune"
                  target="_blank">http://lists.dune-project.org/mailman/listinfo/dune</a><br>
              </div>
            </div>
          </blockquote>
        </div>
        <br>
      </div>
    </blockquote>
    <br>
  </body>
</html>