<div dir="ltr"><div dir="ltr">On Wed, Apr 27, 2022 at 11:31 PM Mike Michell <<a href="mailto:mi.mike1021@gmail.com">mi.mike1021@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">Thanks for the answers. Major reason that not to store those dual cells is to share the same grid file with other codes that are required for coupled physics modeling.<br></div></blockquote><div><br></div><div>Yes, I would create the dual on startup and throw away the original mesh.</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">As a follow-up question, there is struggling to use PetscSection. I attached my example code. <br>The PETSc-provided example, "ex1f90.F90" uses DMPlexCreateSection() to make a PetscSection object. <br>However, by reading .msh & creating a DMPlex, none of label is working with that, except for 'marker' label, which should be provided through "-dm_plex_gmsh_use_marker" option in case gmsh file is used.<br></div></blockquote><div><br></div><div>I have replaced your startup code with the canonical code we use in C examples. It should be more flexible and robust.</div><div>Here is the concise output for your mesh:</div><div><br></div><div>master *:~/Downloads/tmp/Mitchell/Question_3$ ./test -dm_plex_filename new.msh -dm_view <br></div><div>DM Object: 1 MPI processes<br>  type: plex<br>DM_0x84000000_1 in 2 dimensions:<br>  Number of 0-cells per rank: 44<br>  Number of 1-cells per rank: 109<br>  Number of 2-cells per rank: 66<br>Labels:<br>  celltype: 3 strata with value/size (0 (44), 3 (66), 1 (109))<br>  depth: 3 strata with value/size (0 (44), 1 (109), 2 (66))<br>  Face Sets: 4 strata with value/size (5 (5), 7 (5), 6 (5), 8 (5))<br></div><div><br></div><div>The original code was built to take GMsh markers to canonical names because users wanted that. However, if you want your names preserved, you can use</div><div><br></div><div>master *:~/Downloads/tmp/Mitchell/Question_3$ ./test -dm_plex_filename new.msh -dm_view -dm_plex_gmsh_use_regions<br> MPI information ...            0           1<br>DM Object: 1 MPI processes<br>  type: plex<br>DM_0x84000000_1 in 2 dimensions:<br>  Number of 0-cells per rank: 44<br>  Number of 1-cells per rank: 109<br>  Number of 2-cells per rank: 66<br>Labels:<br>  celltype: 3 strata with value/size (0 (44), 3 (66), 1 (109))<br>  depth: 3 strata with value/size (0 (44), 1 (109), 2 (66))<br>  inlet: 1 strata with value/size (5 (5))<br>  upper: 1 strata with value/size (7 (5))<br>  outlet: 1 strata with value/size (6 (5))<br>  lower: 1 strata with value/size (8 (5))<br><br></div><div>There is also code to mark vertices, which some users do not want, that you can turn on using -dm_plex_gmsh_mark_vertices. I use both of these options with PyLith.</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">In my .msh file, there are several boundary conditions, and if I declare labels using them and try to give them to DMGetStratumIS(), the code crashed. It only works with 'marker' label. Any comments?<br>More importantly, trying to use DMPlexCreateSection() is quite painful. </div></blockquote><div><br></div><div>What is not working for you?</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">Therefore, I tried to create PetscSection object through PetscSectionCreate() as attached. But declared vectors are not printed out if I use this way.</div></blockquote><div><br></div><div>If you used the default viewer (PETSC_VIEWER_STDOUT_WORLD), you can see them printed. You are using VTK with a Section that put values on all mesh points. VTK only understands cell fields and vertex fields, so the code does not know how to output this vector. If you define values only on cells or vertices, you will see the output. I believe.</div><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 dir="ltr"> I guess the section object is not properly linked to dm object. Basically, my vector objects (x, y, vel) are not seen from dm viewer & relevant output file. Do you have any recommendations? <br><br>Thanks,<br>Mike<br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">2022년 4월 26일 (화) 오후 6:33, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>>님이 작성:<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, Apr 26, 2022 at 7:27 PM Mike Michell <<a href="mailto:mi.mike1021@gmail.com" target="_blank">mi.mike1021@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">Thanks for the answers. One last question related to designing my code. <br>What I want to do is to build a finite-volume code discretized on "median dual" unstructured mesh. So basically it will use a vertex-centered discretization scheme by looping over "median dual faces", which is not physically existing in my mesh file. <br> <br>I can acheive this by allocating some vectors to compute required geometrical metrics (e.g., face normal vectors, face area, etc. that required to integrate over median dual control volume) as I did before even without PETSc. However, I think this definitely cannot be an optimal way, and I believe there should be much larger benefit that can be obtained, if I well design my DMPlex's data structure (e.g., PetscSection). <br>Is there any example with this median dual control volume case? or Any guideline will be helpful. <br>I was thinking to borrow some functions designed for finite element.<br></div></blockquote><div><br></div><div>I have never done a dual scheme before, so let me ask some questions to better understand the underlying ideas.</div><div><br></div><div>1) Is there a reason not to store the dual directly? That seems like the easiest approach.</div><div><br></div><div>2) You can use Plex to layout auxiliary geometric data. I do it the following way:</div><div><br></div><div>  DMClone(dm, dmGeom);</div><div>  <Create the Section s you want></div><div>  DMSetLocalSection(dmGeom, s);</div><div>  PetscSectionDestroy(&s);</div><div><br></div><div>then</div><div><br></div><div>  DMCreateLocalVector(dmGeom, &lv);</div><div><br></div><div>will give you a Vec with the local data in it that can be addressed by mesh point (global vectors too). Also you would be able</div><div>to communicate this data if the mesh is redistributed, or replicate this data if you overlap cells in parallel.</div><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 dir="ltr">Thanks,<br>Mike<br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">2022년 4월 26일 (화) 오후 5:40, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>>님이 작성:<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, Apr 26, 2022 at 4:48 PM Mike Michell <<a href="mailto:mi.mike1021@gmail.com" target="_blank">mi.mike1021@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">Below two interesting things are found:<br>- If DMPlexCreateGmshFromFile() is used, DMPlexDistribute() is not done by default. I should call DMPlexDistribute() to distribute DMPlex over the procs.<br></div></blockquote><div><br></div><div>Here is the explanation. The default distribution comes from DMSetFromOptions(), which I am guessing was not called after DMPlexCreateGmshFromFile().</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">- If {DMCreate(), DMSetType(), DMSetFromOptions()} in serially used with "-dm_plex_filename ./new.msh" option in command line, DMPlexDistribute() is done by default even though DMPlex object is created by reading the same .msh file.<br><br>Thanks for the modification to the example ex1f90.F90, now it works with 2 procs. <br>But still, if I print my output to "sol.vtk", the file has only a part of whole mesh. However, the file is okay if I print to "sol.vtu". <br>From "sol.vtu" I can see the entire field with rank. Is using .vtu format preferred by petsc?</div></blockquote><div><br></div><div>VTK is generally for debugging, but it should work. I will take a look.</div><div><br></div><div>VTU and HDF5 are the preferred formats.</div><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 class="gmail_quote"><div dir="ltr" class="gmail_attr">2022년 4월 26일 (화) 오후 2:26, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>>님이 작성:<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, Apr 26, 2022 at 9:33 AM Mike Michell <<a href="mailto:mi.mike1021@gmail.com" target="_blank">mi.mike1021@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">Thank you for the answers. <br>For the first question, basically, I cannot run "/dm/impls/plex/ex1f90.F90" example with more than 1 proc. I removed DMPlexDistribute() following your comment and what I tried is: <br><br>- no modification to ex1f90.F90 (as it is)<br>- make "ex1f90"<br>- mpirun -np 2 ./ex1f90<br><br>It gives me "Bad termination of one of ..." for Rank 1. The code runs okay with "mpirun -np 1 ./ex1f90".<br></div></blockquote><div><br></div><div>You are correct. It evidently never worked in parallel, since those checks will only work in serial.</div><div>I have fixed the code, and added a parallel test. I have attached the new file, but it is also in this MR:</div><div><br></div><div>  <a href="https://gitlab.com/petsc/petsc/-/merge_requests/5173" target="_blank">https://gitlab.com/petsc/petsc/-/merge_requests/5173</a></div><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 dir="ltr">Thanks,<br>Mike<br></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">2022년 4월 26일 (화) 오전 3:56, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>>님이 작성:<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 Mon, Apr 25, 2022 at 9:41 PM Mike Michell <<a href="mailto:mi.mike1021@gmail.com" target="_blank">mi.mike1021@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 PETSc developer team,<br><br>I'm trying to learn DMPlex to build a parallel finite volume code in 2D & 3D. More specifically, I want to read a grid from .msh file by Gmsh. <br>For practice, I modified /dm/impls/plex/ex1f90.F90 case to read & distribute my sample 2D grid, which is attached. I have two questions as below:<br><br>(1) First, if I do not use my grid, but use the default box grid built by ex1f90.F90 and if I try "DMPlexDistribute" over mpi processors, the output file (sol.vtk) has only some portion of the entire mesh. How can I print out the entire thing into a single file? Is there any example for parallel output? (Related files attached to "/Question_1/") Paraview gives me some error messages about data size mismatching.<br></div></blockquote><div><br></div><div>For the last release, we made parallel distribution the default. Thus, you do not need to call DMPlexDIstribute() explicitly here. Taking it out, I can run your example.</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">(2) If I create DMPlex object through "DMPlexCreateGmshFromFile", "DMPlexCreateSection" part is crashed. I do not understand why my example code does not work, because the only change was switching from "DMCreate" to "DMPlexCreateGmshFromFile" and providing "new.msh" file. Without the PetscSection object, the code works fine. Any comments about this? (Related files attached to "/Question_2/")<br></div></blockquote><div><br></div><div>If I remove DMPlexDistribute() from this code, it is clear that the problem is with the "marker" label. We do not create this by default from GMsh since we assume people have defined their own labels. You can pass</div><div><br></div><div>  -dm_plex_gmsh_use_marker</div><div><br></div><div>to your code. WHen I do this, you example runs for me.</div><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 dir="ltr">Thanks,<br>Mike<br></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>
</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>
</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>
</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>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><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>