Hi Wayne, it seems that you are using multigrid on a locally refined mesh. Have you tried to use globally refined meshes? This one does not have inner boundaries so that the issue should not occur. Does it?
You might have forgotten to attach the constrains to the transfer operator: https://www.dealii.org/developer/doxygen/deal.II/classMGTransferBlockMatrixFree.html#afd5ace82516e421e37af4572aabf4098. Also you might have forgotten to apply edge constrains. This is how we do it in one of our codes: https://github.com/peterrum/dealii-multigrid/blob/c50581883c0dbe35c83132699e6de40da9b1b255/include/operator.h#L191-L226. Did you set the edge matrices: https://www.dealii.org/developer/doxygen/deal.II/classMultigrid.html#a6d1eff544f75da8c308fdc42c6c0c74d? As proof of concept, you could also try out the global-coarsening multigrid infrastructure (step-75); there, there are no internally boundaries per definition. Peter On Tuesday, 14 February 2023 at 09:25:28 UTC+1 yy.wayne wrote: > Thanks for your reply Martin. > > The reason my Chebyshev preconditioner breaks is the (inverse) diagonal > vector I get from MatrixFree operator contains negative entries. > In my corresponding matrix-based code with exactly same geometry and > boundary condition, the matrix diagonal has no negaive entries. > > I apply a "compute_nonzero_tangential_flux_constraints_on_level" function > similar to the normal_flux one. The geometry of this problem is a shell. > When that BC is applied on outer boundary there's no negative entries, but > the error comes when applied on inner boundary. I'm still debugging to see > where the error comes from. > > I got a quesiton that do FEEvaluation auto deal with constraints? It seems > a matrixfree_operator.vmult solution > equals matrixbased_matrix.vmult + AffineConstraints.distribute. > > Best > Wayne > > 在2023年2月14日星期二 UTC+8 15:33:18<Martin Kronbichler> 写道: > >> To extend over what Wolfgang said, the most likely causes for failure >> with the Chebyshev preconditioner are: >> >> - the vector supplying the matrix diagonal (via the field >> AdditionalData::preconditioner) contains Inf or NaN entries, which in turn >> result from either a ill-formed differential operator or a bug in the code >> computing the diagonal, or >> >> - a coarse level holding only constrained degrees of freedom, which are >> not treated properly so that the SolverCG making an eigenvalue estimate >> needs to work with an unsolvable linear system (right hand side not in the >> image of the matrix). >> >> For the first case, you would check the diagonal on each level and the >> code leading to that, in the second the treatment of constraints. >> >> Best, >> Martin >> On 14.02.23 06:33, 'yy.wayne' via deal.II User Group wrote: >> >> Thank you Wolfgang. I didin't dive much into the Chebyshev preconditioner >> because I'm not familiar with it. >> I assume it's because the CG solver inside Chebyshev preconditioner >> didn't converge or what. >> But as you stated, if NaN should only get from NaN elements and divided >> by zero, then I may find the error. >> I'll follow your debug suggestion. >> >> Best >> Wayne >> >> 在2023年2月14日星期二 UTC+8 12:48:42<Wolfgang Bangerth> 写道: >> >>> >>> > It's strange that the iteration breaks at step 0. If I try right >>> > preconditioning, it breaks at step 1 with nan. Besides, if there is >>> only on >>> > multigrid level, the problem is solved correctly with 1 iteration, so >>> matrices >>> > should be correct. The coarse solution is good enough for smoothing. >>> > >>> > The question is, in which case may a CG or GMRES solver breaks at step >>> 0 with nan? >>> >>> NaN errors, like segmentation faults, are relatively easy to find >>> because you >>> know that every NaN you get is wrong. So you check that the matrix and >>> right >>> hand side and solution vector you give to the solver don't have any >>> NaNs. If >>> they don't, and you get NaNs out, then the problem likely lies with the >>> preconditioner, so you single step through the solver around the place >>> where >>> the preconditioner is applied and check in which operation the NaN >>> appears. >>> >>> In your case, you seem to have already found out that it's the Chebyshev >>> preconditioner. The question is why. Is it because the matrix the >>> preconditioner uses has NaNs in it? Is it because the preconditioner >>> makes an >>> assumption that it can divide by certain matrix entries, but they are >>> zero? >>> Single stepping in a debugger will give you the answer :-) >>> >>> Best >>> W. >>> >>> -- >>> ------------------------------------------------------------------------ >>> Wolfgang Bangerth email: bang...@colostate.edu >>> www: http://www.math.colostate.edu/~bangerth/ >>> >>> >>> -- >> The deal.II project is located at http://www.dealii.org/ >> For mailing list/forum options, see >> https://groups.google.com/d/forum/dealii?hl=en >> --- >> You received this message because you are subscribed to the Google Groups >> "deal.II User Group" group. >> To unsubscribe from this group and stop receiving emails from it, send an >> email to dealii+un...@googlegroups.com. >> To view this discussion on the web visit >> https://groups.google.com/d/msgid/dealii/a0326f84-a7b7-49ec-97b7-db1579574e87n%40googlegroups.com >> >> <https://groups.google.com/d/msgid/dealii/a0326f84-a7b7-49ec-97b7-db1579574e87n%40googlegroups.com?utm_medium=email&utm_source=footer> >> . >> >> -- The deal.II project is located at http://www.dealii.org/ For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en --- You received this message because you are subscribed to the Google Groups "deal.II User Group" group. To unsubscribe from this group and stop receiving emails from it, send an email to dealii+unsubscr...@googlegroups.com. To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/acddfe4d-29cb-4d1b-8c8a-509657f234een%40googlegroups.com.