<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Jul 15, 2016 at 10:48 PM, Barry Smith <span dir="ltr"><<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><br>
> On Jul 14, 2016, at 12:21 PM, domenico lahaye <<a href="mailto:domenico_lahaye@yahoo.com">domenico_lahaye@yahoo.com</a>> wrote:<br>
><br>
> Dear PETSc team,<br>
><br>
> 1) I am looking into ks/examples/tutorials/ex42.c<br>
<br>
This example is really written as only a one level solver, making it work with geometric multigrid is not clean<br>
<br>
> I am still new to the DMDA structure<br>
> and likely not giving it as much time as it deserves. However, I do not see immediately<br>
> what function is responsible for calling PCMGSetSmoother and PCMGSetResidual.<br>
><br>
> I tried to call PCMGGetCoarseSolve(pc, &kcpc) and subsequently<br>
> KSPGetOperators (kspc, ... ) to check how the coarse grid operator is defined<br>
> after calling DMCoarsenHierarchy, but that failed.<br>
><br>
> I am solving Helmholtz with shifted Laplace, and managed to exploit DMDA to perform<br>
> a multigrid solve on the preconditioner. In a next stage I want to implement the deflation<br>
> using DMDA as well.<br>
<br>
You should look at ex25.c in the same directory. Here<br>
<br>
ierr = KSPSetDM(ksp,da);CHKERRQ(ierr);<br>
ierr = KSPSetComputeRHS(ksp,ComputeRHS,&user);CHKERRQ(ierr);<br>
ierr = KSPSetComputeOperators(ksp,ComputeMatrix,&user);CHKERRQ(ierr);<br>
<br>
make it straight forward to work with multigrid. The KSP object can mange the hierarchy of grids since it is provided with the DM<br>
and the ComputeRHS and ComputeMatrix provide a way for the multigrid preconditioner to automatically generate the needed matrix on each level without you having to manage it yourself. For example the rule in the makefile<br>
<br>
runex25:<br>
-@${MPIEXEC} -n 1 ./ex25 -pc_type mg -ksp_type fgmres -da_refine 2 -ksp_monitor_short -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -ksp_view -pc_mg_type full > ex25_1.tmp 2>&1; \<br>
if (${DIFF} output/ex25_1.out ex25_1.tmp) then true; \<br>
else printf "${PWD}\nPossible problem with ex25_1, diffs above\n=========================================\n"; fi; \<br>
${RM} -f ex25_1.tmp<br>
<br>
shows how to run with two levels. etc.<br>
<br>
<br>
><br>
> 2) On <a href="http://www.mcs.anl.gov/petsc/documentation/referencing.html" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/documentation/referencing.html</a> I see<br>
><br>
> @Misc{petsc-web-page,<br>
> author = {Satish Balay and Shrirang Abhyankar and Mark~F. Adams and Jed Brown and Peter Brune<br>
> and Kris Buschelman and Lisandro Dalcin and Victor Eijkhout and William~D. Gropp<br>
> and Dinesh Kaushik and Matthew~G. Knepley<br>
> and Lois Curfman McInnes and Karl Rupp and Barry~F. Smith<br>
> and Stefano Zampini and Hong Zhang and Hong Zhang},<br>
> title = {{PETS}c {W}eb page},<br>
> url = {<a href="http://www.mcs.anl.gov/petsc" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc</a>},<br>
> howpublished = {\url{<a href="http://www.mcs.anl.gov/petsc}" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc}</a>},<br>
> year = {2016}<br>
> }<br>
><br>
><br>
><br>
> Is the last author mentioned twice intentionally?<br></blockquote><div><br></div><div>That is actually two different people with the same name.</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
> 3) On <a href="http://www.mcs.anl.gov/petsc/publications/petscapps-bib.html#OpenFOAM%202.2.1" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/publications/petscapps-bib.html#OpenFOAM%202.2.1</a> I see<br>
><br>
> @misc{OpenFOAM<br>
> ,<br>
><br>
><br>
> title = "OpenFOAM",<br>
><br>
> howpublished = "\url{<a href="http://www.openfoam.com" rel="noreferrer" target="_blank">http://www.openfoam.com</a>}",<br>
><br>
> url = {<a href="http://www.openfoam.com" rel="noreferrer" target="_blank">http://www.openfoam.com</a>},<br>
><br>
> note = "OpenFOAM is a free, open source CFD software package. It allows PETSc linear algebra and solvers to be used underneath.",<br>
><br>
> key = "OpenFOAM 2.2.1"<br>
><br>
> }<br>
><br>
><br>
> Do you have more information on the use of PETSc within OpenFoam?<br></blockquote><div><br></div><div>They only use solvers, and not the DM stuff as far as I know.</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
> 4) @matt in response to a question he raised in Vienna<br>
><br>
> MIPSE is a BEM solver. Details are on:<br>
> <a href="http://www.g2elab.grenoble-inp.fr/plateforms/mipse-modeling-of-interconnected-power-systems-632862.kjsp?RH=G2ELAB_R-MAGE" rel="noreferrer" target="_blank">http://www.g2elab.grenoble-inp.fr/plateforms/mipse-modeling-of-interconnected-power-systems-632862.kjsp?RH=G2ELAB_R-MAGE</a></blockquote><div><br></div><div>From what I can tell, the code is not open source. Is that right?</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><br>
> Cheers, Domenico Lahaye.<br>
><br>
<br>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature" data-smartmail="gmail_signature">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></div>