<div dir="ltr"><div dir="ltr">On Wed, Nov 10, 2021 at 8:26 AM Bhargav Subramanya <<a href="mailto:bhargav.subramanya@kaust.edu.sa">bhargav.subramanya@kaust.edu.sa</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 Matt,<div><br><div>Thanks a lot for the reply. I am now able to generate the prismatic mesh properly. </div></div></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 dir="ltr"><div><div>In the case of mpiexec -n 2 ./meshtest -dm_plex_shape sphere -dm_refine_pre 3 -dm_distribute -dm_refine 2 and DMExtrude(), where I am ordering the extruded cells on the layers first; Is the mesh extruded in parallel for this case? </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 dir="ltr"><div><div>The first time I distributed the mesh, it was only the spherical surface mesh. After dmplex extrude, the mesh changes. Do I need to redistribute the mesh again?</div></div></div></blockquote><div><br></div><div>Probably not. The balance factor will not change, meaning the relative sizes of the partitions. However, if you have a lot of layers, maybe you want to. It is easy,</div><div>just call DMDistribute() again.</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"><div><div>Thanks,</div><div>Bhargav</div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Nov 10, 2021 at 12:56 AM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:<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, Nov 9, 2021 at 9:54 AM Bhargav Subramanya <<a href="mailto:bhargav.subramanya@kaust.edu.sa" target="_blank">bhargav.subramanya@kaust.edu.sa</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 want to generate a prismatic mesh from a base spherical surface mesh. I first refine the base icosahedron sequential mesh to some extent, distribute the mesh, continue the refinement in parallel, project the mesh to a unit sphere, and then finally extrude the mesh. The following is the code I am using. As the extrusion is performed in parallel, the extruded mesh seems to be broken as shown in the attached figure. May I know how to get an intact extruded mesh in parallel? Also, is it possible to make mesh refinement respect the spherical surface geometry, without having to project it using a function as shown below?</div></div></blockquote><div><br></div><div>I think you can do most of what you want from the command line. Here is a simple code:</div><div><br></div><div>#include <petsc.h><br><br>int main(int argc, char **argv)<br>{<br>  DM             dm;<br>  PetscErrorCode ierr;<br><br>  ierr = PetscInitialize(&argc, &argv, NULL, NULL);if (ierr) return ierr;<br>  ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);<br>  ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);<br>  ierr = DMSetFromOptions(dm);CHKERRQ(ierr);<br>  ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);<br>  ierr = DMDestroy(&dm);CHKERRQ(ierr);<br>  ierr = PetscFinalize();<br>  return ierr;<br>}<br></div><div><br></div><div>that I run using</div><div><br></div><div>  mpiexec -n 2 ./meshtest -dm_plex_shape sphere -dm_refine_pre 3 -dm_extrude 5 -dm_plex_transform_extrude_thickness 0.1 -dm_distribute -dm_view hdf5:mesh.h5</div><div><br></div><div>and get the attached picture. I am refining before distribution only here, since the after distribution refinement would refine the extrude thing, which is not what you want.</div><div>If you want refinement after distribution, you could use this code to do</div><div><br></div><div>  mpiexec -n 2 ./meshtest -dm_plex_shape sphere -dm_refine_pre 3 -dm_distribute -dm_refine 2</div><div><br></div><div>and then call DMExtrude() by hand after that. Alternatively, I could make a pre-extrusion refinement option.</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"><div><div><div>Thanks,</div><div>Bhargav</div><div><br></div><div>/* refine the sequential mesh first */<br>                  for (int r = 0; r < numRefinePre; ++r) {<br>                          DM rdm = NULL;<br>                                DMRefine(dm, MPI_COMM_WORLD, &rdm);<br>                               DMDestroy(&dm);<br>                           dm = rdm;<br>                       }<br><br>                         /* project to a unit sphere */<br>                        ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);<br><br>                      /* create and distribute DM */<br>                        ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);<br>                      if (dmDist) {<br>                                DMDestroy(&dm);<br>                           dm   = dmDist;<br>                        }<br><br>                   /* refine the mesh in parallel */<br>                     for (int r = 0; r < numRefine; ++r) {<br>                             DM rdm = NULL;<br>                                DMRefine(dm, MPI_COMM_WORLD, &rdm);<br>                               DMDestroy(&dm);<br>                           dm = rdm;<br>                       }<br><br>                         /* project to a unit sphere */<br>                        ierr = ProjectToUnitSphere(dm); CHKERRQ(ierr);<br><br>                      ierr = DMPlexExtrude(dm, numLayers-1, radialThickness, PETSC_FALSE, NULL, PETSC_TRUE, &dm3D);CHKERRQ(ierr);<br><br>                     if (dm3D) {<br>                          DMDestroy(&dm);<br>                           dm = dm3D;<br>                     }<br></div><div><br></div></div></div></div>

<br>
<div><hr></div><font face="Arial" size="1">This message and its contents, including attachments are intended solely for the original recipient. If you are not the intended recipient or have received this message in error, please notify me immediately and delete this message from your computer system. Any unauthorized use or distribution is prohibited. Please consider the environment before printing this email.</font></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="https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!Nmw4Hv0!lx1-mURFRmaF1Jo4rDYcZmLsnPsG8EO1HFWoGa2ukYZRt1fF4JLOpS3gHpaUHzL6_KersFjJlmA$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>
</blockquote></div>

<br>
<div><hr></div><font face="Arial" size="1">This message and its contents, including attachments are intended solely for the original recipient. If you are not the intended recipient or have received this message in error, please notify me immediately and delete this message from your computer system. Any unauthorized use or distribution is prohibited. Please consider the environment before printing this email.</font></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>