<div dir="ltr"><div dir="ltr">On Mon, Mar 11, 2019 at 10:42 AM Pietro Benedusi via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</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 style="word-wrap:break-word">
Dear Petsc team,
<div><br>
</div>
<div>I have a question about the setting up of a multigrid solver.</div>
<div><br>
</div>
<div>I would like yo use a PCG smoother, preconditioned with a mass matrix, just on the fine level.</div>
<div>But when add the line for preconditioning the CG with the mass matrix my MG diverges.</div>
<div><br>
</div>
<div>I have implemented the same solver in MATLAB and it converges fine. Also the operators in PETSc are the same and the PCG applied directly on the problem (without MG) works the same in both PETSC and MATLAB. </div>
<div><br>
</div>
<div>This is what I do in PETSC for 2 levels:</div></div></blockquote><div><br></div><div>Send the output of -ksp_view so we can see what the exact solver is. It will not tell whether the matrices you set are the correct ones, but at least we can see the solver.</div><div><br></div><div>Also, did you look at the residual decrease? Is  the smoother working on every level but the finest? Is the coarse grid correction effective?</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 style="word-wrap:break-word">
<div><font face="Menlo"><br>
</font></div>
<div>
<div><font face="Menlo">    KSP space_solver;       </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">    ierr = KSPCreate(PETSC_COMM_WORLD,&space_solver);CHKERRQ(ierr);</font></div>
<div><font face="Menlo">    ierr = KSPSetTolerances(space_solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);           </font></div>
<div><font face="Menlo">    ierr = KSPSetOperators(space_solver, K, K);CHKERRQ(ierr); </font></div>
<div><font face="Menlo">    ierr = KSPSetNormType(space_solver,  KSP_NORM_UNPRECONDITIONED );CHKERRQ(ierr);  </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">    ierr = KSPSetType(space_solver,KSPRICHARDSON);CHKERRQ(ierr);    </font></div>
<div><font face="Menlo">    ierr = KSPSetFromOptions(space_solver);CHKERRQ(ierr);         </font></div>
<div><font face="Menlo">    ierr = KSPSetUp(space_solver);CHKERRQ(ierr);        </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">    PC pcmg;</font></div>
<div><font face="Menlo">    ierr = KSPGetPC(space_solver, &pcmg);</font></div>
<div><font face="Menlo">    ierr = PCSetType(pcmg, PCMG);    </font></div>
<div><font face="Menlo">    ierr = PCMGSetLevels(pcmg,levels, NULL);CHKERRQ(ierr);</font></div>
<div><font face="Menlo">    ierr = PCMGSetGalerkin(pcmg,PC_MG_GALERKIN_BOTH);CHKERRQ(ierr);</font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">    // smoothers</font></div>
<div><font face="Menlo">    for (int i = 1; i < levels; ++i)  </font></div>
<div><font face="Menlo">    {</font></div>
<div><font face="Menlo">        KSP smoother;</font></div>
<div><font face="Menlo">        ierr = PCMGGetSmoother(pcmg, i, &smoother);CHKERRQ(ierr);</font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">            ierr = KSPSetType(smoother, KSPCG);CHKERRQ(ierr); </font></div>
<div><font face="Menlo">            ierr = KSPSetOperators(smoother, K, M);CHKERRQ(ierr);             </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">            // ierr = KSPSetUp(smoother);CHKERRQ(ierr);              </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">            ierr = KSPSetTolerances(smoother,1e-12,PETSC_DEFAULT,PETSC_DEFAULT,s_p);CHKERRQ(ierr);               </font></div>
<div><font face="Menlo">            ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);  </font></div>
<div><font face="Menlo"><br>
</font></div>
<div><font face="Menlo">            PC sm;</font></div>
<div><font face="Menlo">            ierr = KSPGetPC(smoother, &sm);CHKERRQ(ierr);            </font></div>
<div><font face="Menlo">            ierr = PCSetType(sm,   PCLU);CHKERRQ(ierr);       </font></div>
<div><font face="Menlo">                        </font></div>
<div><font face="Menlo">        ierr = PCMGSetInterpolation(pcmg, i, interpolation_operators[i-1]);CHKERRQ(ierr);                   </font></div>
<div><font face="Menlo">    }</font></div>
<div><br>
</div>
<div><br>
</div>
<div><br>
</div>
<div>I think there is a problem with the PETSc syntax, because I checked everything else and it is fine.</div>
<div><br>
</div>
<div>Do you any ideas?</div>
<div><br>
</div>
<div>Thank you very much!</div>
<div><br>
</div>
<div>Best,</div>
<div>Pietro </div>
<div><br class="gmail-m_851211825493033580webkit-block-placeholder">
</div>
<div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
~~~~~~~~~~~~</div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
Pietro Benedusi</div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
<br>
</div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
Numerical Simulation in Science,<br>
Medicine and Engineering research group</div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
ICS, Institute of Computational Science<br>
</div>
<div style="color:rgb(0,0,0);letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;word-wrap:break-word">
USI, Università della Svizzera Italiana<br>
Via Giuseppe Buffi, 13<br>
CH - 6900 Lugano<br>
<div><a href="mailto:benedp@usi.ch" target="_blank">benedp@usi.ch</a><br>
<br>
</div>
</div>
</div>
</div>
</div>
<br>
</div>
</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>