<div dir="ltr"><div>Yup -- I realized that I had '-pc_type mg' in the script I was using to build and run as I developed. I guess that was causing the KSP to coarsen its DM, which made me think I had to force the density DM to be consistent. Still, refactoring my original grid DM into one for Vlasov components and one for the potential was useful.</div><div><br></div><div>Thanks for the help!<br></div><div><div><div><div><div dir="ltr" class="gmail_signature" data-smartmail="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div><br></div><div>--Matt</div><div dir="ltr">==========================<br></div><div dir="ltr">Matthew Young, PhD (he/him)<div><div><div><div>Research Scientist II</div><div>Space Science Center<br></div></div></div><div>University of New Hampshire</div><div><a href="mailto:Matthew.Young@unh.edu" target="_blank">Matthew.Young@unh.edu</a><br></div></div><div><div dir="ltr"><div dir="ltr">==========================<br></div></div></div></div></div></div></div></div></div></div><br></div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Mon, May 1, 2023 at 10:51 PM Dave May <<a href="mailto:dave.mayhem23@gmail.com">dave.mayhem23@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><br></div><div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Mon 1. May 2023 at 18:57, Matthew Young <<a href="mailto:myoung.space.science@gmail.com" target="_blank">myoung.space.science@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>Thanks for the suggestion to keep DMs separate, and for pointing me toward that example. I now have a DM for the particle quantities (i.e., density and flux) and another for the potential. I'm hoping to use KSPSetComputeOperators with PCGAMG, so I packed the density DM into the application context and set the potential DM on the KSP, but I'm not sure how to communicate changes in the KSP DM (e.g., coarsening) to the density DM inside my operator function.</div></div></blockquote><div dir="auto"><br></div><div dir="auto">I don’t think you need to.</div><div dir="auto"><br></div><div dir="auto">GAMG only requires the fine grid operator - this will be the matrix assembled from KSPSetComputeOperators. Hence density DM and potential DM fields only need to be managed by you on the finest level. </div><div dir="auto"><br></div><div dir="auto">However, if you wanted to use PCMG with rediscretized operators on every level, then you would need the density DM field defined on each level of your geometric multigrid hierarchy. This could be done (possibly less than ideally) by calling DMCreateInterpolation() and then using the Mat to interpolate the density from the  finest level to next coarsest level (and so on).</div><div dir="auto"><br></div><div dir="auto">Thanks,</div><div dir="auto">Dave </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 dir="auto"><br></div></div><div dir="ltr"><div><div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div><br></div><div>--Matt</div><div dir="ltr">==========================<br></div><div dir="ltr">Matthew Young, PhD (he/him)<div><div><div><div>Research Scientist II</div><div>Space Science Center<br></div></div></div><div>University of New Hampshire</div><div><a href="mailto:Matthew.Young@unh.edu" target="_blank">Matthew.Young@unh.edu</a><br></div></div><div><div dir="ltr"><div dir="ltr">==========================<br></div></div></div></div></div></div></div></div></div></div><br></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sun, Apr 30, 2023 at 1:52 PM 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 Sun, Apr 30, 2023 at 1:12 PM Matthew Young <<a href="mailto:myoung.space.science@gmail.com" target="_blank">myoung.space.science@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"><div>Hi all,</div><div><br></div><div>I am developing a particle-in-cell code that models ions as particles and electrons as an inertialess fluid. I use a PIC DMSWARM for the ions, which I gather into density and flux before solving a linear system for the electrostatic potential (phi). I currently have one DMDA with 5 degrees of freedom -- one each for density, 3 flux components, and phi.</div><div><br></div><div>When setting up the linear system to solve for phi, I've been following examples like KSP ex34.c and ex42.c when writing the KSP operator and RHS functions but I'm not sure I have the right approach, since 4 of the DOFs are known and 1 is unknown.</div><div><br></div><div>I saw <a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2016-November/031031.html" target="_blank">this thread</a> that recommended using DMDAGetReducedDMDA, which I gather has been deprecated in favor of DMDACreateCompatibleDMDA. Is that a good approach for managing a regular grid with known and unknown quantities on each node? Could a composite DM be useful? Has anyone else worked on a problem like this?</div></div></blockquote><div><br></div><div>I recommend making a different DM for each kind of solve you want. DMDACreateCompatibleDMDA() should be the implementation of DMClone(), but we have yet to harmonize all things for all DMs. I would create one DM for your Vlasov components and one for the Poisson.</div><div>We follow this strategy in our Vlasov-Poisson test for Landau damping: <a href="https://gitlab.com/petsc/petsc/-/blob/main/src/dm/impls/swarm/tests/ex9.c" target="_blank">https://gitlab.com/petsc/petsc/-/blob/main/src/dm/impls/swarm/tests/ex9.c</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"><div><div><div dir="ltr"><div dir="ltr"><div><div dir="ltr"><div><div>--Matt</div><div dir="ltr">==========================<br></div><div dir="ltr">Matthew Young, PhD (he/him)<div><div><div><div>Research Scientist II</div><div>Space Science Center<br></div></div></div><div>University of New Hampshire</div><div><a href="mailto:Matthew.Young@unh.edu" target="_blank">Matthew.Young@unh.edu</a><br></div></div><div><div dir="ltr"><div dir="ltr">==========================<br></div></div></div></div></div></div></div></div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div><span>-- </span><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></div>
</blockquote></div>