<div dir="ltr"><div dir="ltr">Hello Matt, <div><br></div></div><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le mar. 18 mai 2021 à 20:31, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Tue, May 18, 2021 at 12:16 PM Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="auto">Hi Matthew,</div><div dir="auto"><br></div><div dir="auto">Thank you very much for your quick answer, as always ! </div></blockquote><div><br></div><div>Cool.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le mar. 18 mai 2021 à 17:46, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> a écrit :<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div dir="ltr">On Tue, May 18, 2021 at 5:19 AM Thibault Bridel-Bertomeu <<a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr">Dear all,<div><br></div><div>I am playing around with creating DMPlex from a periodic Gmsh mesh (still in the same finite volume code that solves the Euler equations) because I want to run some classical periodic test cases from the literature (like convecting a vortex and so on).</div><div>I got from the mailing list and the examples that I gotta use -dm_plex_gmsh_periodic for the reader and the DM creator to do their job properly and read the $Periodic section from the Gmsh mesh. That works fine, I can write the DM in HDF5 format afterwards and the periodicity is clearly there (i can see in Paraview that opposite sides are "attached" to each other).</div><div>Now, in the context of the finite volume solver, I wrote my one right-hand side term computing routine, and I rely heavily on routines such as DMPlexGetFaceFields. Those routines, in turn rely on the ghost cells and ... I am having an issue with DMPlexConstructGhostCells. When I use -dm_plex_gmsh_periodic, DMPlexCreateFromFile works fine, DMSetBasicAdjacency(dm, true, false) works fine, DMPlexDistribute works fine, but when I get to DMPlexConstructGhostCells, I get</div><div><br></div><div> "<span style="font-variant-ligatures:no-common-ligatures;font-family:Menlo;font-size:11px;background-color:rgb(252,238,207);color:rgb(0,0,0)">DM has boundary face 2048 with 2 support cells"</span></div>





<div><br clear="all"><div><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div><div><div>I understand the DM begin periodic, even the "boundary" faces have two neighbors because of the topological folding, so ... yeah, the ghost cells cannot really exist. But then how can I proceed ? Because I need that "ghost" label for all the routines I call from the RHS computing routine ...</div></div></div></div></div></div></div></div></div></div></blockquote><div><br></div><div>Thanks for the example. I have run this myself.  It looks to me like the mesh you are providing. Has no boundary. It is doubly periodic (a torus). When I tell ConstructGhostCells() to</div><div>create a label for the boundary, rather than using "Face Sets", it creates an empty label because no faces have only 1 neighbor.</div></div></div></blockquote><div dir="auto"><br></div><div dir="auto">Yes indeed, up is down and right is left for this mesh, the idea being of emulating an infinite medium. </div><div dir="auto">So you manage to get through ConstructGhostCells without a segfault ? What did you specify instead of PETSC_NULL_CHARACTER ?</div></div></div></blockquote><div><br></div><div>The NULL tells it to use the label "Face Sets", but that label coming from Gmsh marks edges that are not actually boundaries. I asked PETSc</div><div>to create a boundary label named "marker" using <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexMarkBoundaryFaces.html" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexMarkBoundaryFaces.html</a></div><div>It was empty, but everything goes through fine.</div></div></div></blockquote><div><br></div><div>OK thanks ! I had to write the F90 wrapper for DMLabelCreate but aside from that, everything indeed goes through fine with the ConstructGhostCells.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div class="gmail_quote"><div dir="auto">Also, does it still create a « ghost » label even if it’s empty ?</div></div></div></blockquote><div><br></div><div>Yes.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div class="gmail_quote"><div dir="auto"> Do you think that will hold for later uses in the GetFaceFields and similar routines ? </div></div></div></blockquote><div><br></div><div>Yes.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div class="gmail_quote"><div dir="auto">I read through the source and it seems it won’t like it much if there are no ghost at all, but since I did not pass the ConstructGhostCell step I’m not sure. </div></div></div></blockquote><div><br></div><div>It should work fine, you just will not need boundary conditions.</div></div></div></blockquote><div><br></div><div>Actually, there is an issue with the corner cells. My mesh is periodic in both directions and it turns out the corner cells accumulate way too much numerical flux and the computation crashes after 4 iterations or so - it happens both with quads and tris. Whether I use my version of the RHS computing routine or the DMPlexTSComputeRHSFunctionFVM leads to the same crash.</div><div>I'll keep investigating, but if you have any idea what could go wrong ...</div><div><br></div><div>Thanks !!</div><div><br></div><div>Thibault</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><div class="gmail_quote"><div dir="auto">Thanks ! </div><div dir="auto"><br></div><div dir="auto">Thibault</div><div dir="auto"><br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><div dir="auto"></div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div></div></div><div dir="ltr"><div class="gmail_quote"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div><div><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div><div><div>I attach to this email the periodic Gmsh mesh I am playing with, and the routine from my code where I create the DM and so on:</div><div><br></div><div><div style="font-family:Menlo,Monaco,"Courier New",monospace;font-size:12px;line-height:18px;white-space:pre-wrap;background-color:rgb(30,30,30);color:rgb(212,212,212)"><div style="font-family:Menlo,Monaco,"Courier New",monospace">        <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(86,156,214)">subroutine</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">initmesh</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            PetscErrorCode :: ierr</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            DM             :: dmDist, dmGhost</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(86,156,214)">integer</span>        :: overlap</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            PetscViewer    :: vtkViewer</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscLogEventBegin</span>(timer_initmesh, ierr); CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Number of neighbours taken into account in MP communications(1 - Order 1; 2 - Order 2)</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            overlap = <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(181,206,168)">1</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscPrintf</span>(PETSC_COMM_WORLD, <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(206,145,120)">"Initializing mesh...\n"</span>, ierr)          ; CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Force DMPlex to use gmsh marker</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! call PetscOptionsSetValue(PETSC_NULL_OPTIONS, "-dm_plex_gmsh_use_marker", "true", ierr); CHKERRA(ierr)</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Read mesh from file name 'meshname'</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMPlexCreateFromFile</span>(PETSC_COMM_WORLD, meshname, PETSC_TRUE, dm, ierr); CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Distribute on processors</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Start with connectivity</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMSetBasicAdjacency</span>(dm, PETSC_TRUE, PETSC_FALSE, ierr)                 ; CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Distribute on processors</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMPlexDistribute</span>(dm, overlap, PETSC_NULL_SF, dmDist, ierr)             ; CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Security check</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">if</span> (dmDist <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(86,156,214)">/=</span> PETSC_NULL_DM) <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">then</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Destroy previous dm</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMDestroy</span>(dm, ierr)                                                ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Replace with dmDist</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                dm = dmDist</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">end if</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Finalize setup of the object</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMSetFromOptions</span>(dm, ierr)                                             ; CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Boundary condition with ghost cells</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMPlexConstructGhostCells</span>(dm, PETSC_NULL_CHARACTER, PETSC_NULL_INTEGER, dmGhost, ierr); CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Security check</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">if</span> (dmGhost <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(86,156,214)">/=</span> PETSC_NULL_DM) <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">then</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Destroy previous dm</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMDestroy</span>(dm, ierr)                                                ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Replace with dmGhost</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                dm = dmGhost</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">end if</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">if</span> (debug) <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">then</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! Show in terminal</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscPrintf</span>(PETSC_COMM_WORLD, <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(206,145,120)">":: [DEBUG] Visualizing DM in console ::\n"</span>, ierr); CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMView</span>(dm, PETSC_VIEWER_STDOUT_WORLD, ierr)                        ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(106,153,85)">! VTK viewer</span></div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscViewerCreate</span>(PETSC_COMM_WORLD, vtkViewer, ierr)               ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscViewerSetType</span>(vtkViewer, PETSCVIEWERHDF5, ierr)               ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscViewerFileSetMode</span>(vtkViewer, FILE_MODE_WRITE, ierr)           ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscViewerFileSetName</span>(vtkViewer, <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(206,145,120)">"debug_initmesh.h5"</span>, ierr)       ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">DMView</span>(dm, vtkViewer, ierr)                                        ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">                <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscViewerDestroy</span>(vtkViewer, ierr)                                ; CHKERRA(ierr)</div><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">end if</span></div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscPrintf</span>(PETSC_COMM_WORLD, <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(206,145,120)">"Done !\n"</span>, ierr)                        ; CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">            <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(197,134,192)">call</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">PetscLogEventEnd</span>(timer_initmesh, ierr); CHKERRA(ierr)</div><br><div style="font-family:Menlo,Monaco,"Courier New",monospace">        <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(86,156,214)">end subroutine</span> <span style="font-family:Menlo,Monaco,"Courier New",monospace;color:rgb(220,220,170)">initmesh</span></div></div></div><div><br></div><div>Thank you very much for your help !</div><div><br></div><div>Best regards, </div><div><br></div><div>Thibault Bridel-Bertomeu</div></div></div></div></div></div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div dir="ltr"><div><div><div><div><div>Thibault Bridel-Bertomeu<br>—<br></div></div></div></div>Eng, MSc, PhD</div><div>Research Engineer</div><div>CEA/CESTA</div><div>33114 LE BARP</div><div>Tel.: (+33)557046924</div><div>Mob.: (+33)611025322<br></div><div>Mail: <a href="mailto:thibault.bridelbertomeu@gmail.com" target="_blank">thibault.bridelbertomeu@gmail.com</a><br></div></div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div></div>