[petsc-users] [petsc-maint] Poisson with Neumann BC on a stretched grid

Barry Smith bsmith at petsc.dev
Thu Oct 8 21:31:12 CDT 2020



> On Oct 8, 2020, at 2:50 AM, Victoria Hamtiaux <victoria.hamtiaux at uclouvain.be> wrote:
> 
> I'm sorry but I'm a bit confused. 
> 
> 
> 
> First, the fact that the matrix is not symmetric is ok?
> 
> 
> 
> Secondly, I guess that coding the semi-coarsening would be the "best" solution, but isn't there any "easier" solution to solve that linear system (Poisson equation with pure Neumann BC on a stretched grid (and in parallel))? 
> 
> 
> 
> Also, is it normal that using the PCLU preconditionner (and solving on 1 processor) with a stretched grid,  the solution of the linear solver is bad. Is it possible that the PCLU preconditionner also has problems with stretched grids? (Because again, with a uniform grid, it works fine).
> 

  No this is not normal. I would expect PCLU to be very robust for stretched grids. 

  With finite differences I think a stretched grid results in a non symmetric matrix because the differencing (to the right and left for example) uses a different h. I think this is normal, with traditional finite elements I think it will remain symmetric. Also with finite differencing the order of the accuracy of the discretization falls because the terms in the Taylor series that cancel with a non-stretched grid do not cancel with a stretched grid.

   Is it possible that something is wrong with generation of the matrix with the stretched grid?  You could try a convergence study with a slightly stretched grid using PCLU to see if it seems to be working correctly. Only when the numerics are right and you want to run a large problem where PCLU is too slow you would switch to multgrid with semi-coarsening. I think it is pretty easy for a structured grid and we can show you how and maybe get a nice example for PETSc out of it.

  Barry




> 
> Sorry for asking so much questions, 
> 
> 
> 
> Thanks again for your help, 
> 
> 
> 
> Best regards, 
> 
> 
> 
> Victoria
> 
> 
> 
> 
> 
> On 7/10/20 17:59, Barry Smith wrote:
>> 
>> 
>>> On Oct 7, 2020, at 8:11 AM, Mark Adams <mfadams at lbl.gov <mailto:mfadams at lbl.gov>> wrote:
>>> 
>>> 
>>> On Wed, Oct 7, 2020 at 8:27 AM Victoria Hamtiaux <victoria.hamtiaux at uclouvain.be <mailto:victoria.hamtiaux at uclouvain.be>> wrote:
>>> Thanks for all the answers,
>>> 
>>> 
>>> 
>>> How can I do the "semi-coarsening"? I don't really know how those preconditionners work so I don't how how to change them or so..
>>> 
>>> 
>>> 
>>> You have to write custom code to do semi-coarsening. PETSc does not provide that and you would not want to do it yourself, most likely.
>> 
>>   We do not provide it directly but if you are using PCMG and DMDA it is relatively straight-forward. You create a coarse DM and then refine it but each refinement you only do in the directions you want set each time with DMDASetRefinementFactor(). Once you have the collections of refined DM's you provide them to PCMG.
>> 
>>   Barry
>> 
>>>  
>>> 
>>> I have a question because you both seem to say that my matrix is supposed to be symmetric which is not the case. \
>>> 
>>> You said "my matrix is symmetric." 
>>> 
>>> Then you said " I suspect that by stretching the grid, my matrix is not symmetric anymore and that it might cause a problem."
>>> 
>>> We are saying that by stretchin the grid the matrix is still symmetric even if the grid has lost a symmetry. I don't know of a mechanism for stretching the grid to make the matrix asymmetric. So we are suggesting that you verify your suspicion that the matrix is symmetric.
>>> 
>>> And in fact, I don't get how it can be symmetric. Because you will have something close to symmetric. For example when you are at the center of your domain it will be symmetric, but when your at a point at the boundaries I don't get how you can be symmetric, you won't have something at the left and the right of your main diagonal... (I don't know if my explanations are understandable)
>>> 
>>> You can make a discretization that is not symmetric because of boundary conditions but I assume that is not the case because you said your matrix is symmetric.
>>>  
>>> Best regards, 
>>> 
>>> 
>>> 
>>> Victoria
>>> 
>>> 
>>> 
>>> 
>>> 
>>> On 7/10/20 14:20, Mark Adams wrote:
>>>> GMG (geometric MG) is stronger as Matt said, but it is affected by stretched grids in a similar way. A way to fix this in GMG is semi-coarsening, which AMG _can_ do automatically.
>>>> 
>>>> On Wed, Oct 7, 2020 at 8:02 AM Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>>> On Wed, Oct 7, 2020 at 7:07 AM Victoria Hamtiaux <victoria.hamtiaux at uclouvain.be <mailto:victoria.hamtiaux at uclouvain.be>> wrote:
>>>> Hello Matt, 
>>>> 
>>>> 
>>>> 
>>>> I just checked the symmetry of my matrix and it is not symmetric. But it is not symmetric either when I use a uniform grid.
>>>> 
>>>> The domain is 3D and I'm using finite differences, so I guess it is normal that at multiple places (when I deal with points near the boundaries), the matrix is not symmetric.
>>>> 
>>>> So I was wrong, the problem doesn't come from the fact that the matrix is not symmetric. I don't know where it comes from, but when I switch from uniform to stretched grid, the solver stops working properly. Could it be from the preconditionner of the solver that I use?
>>>> 
>>>> Do you have any other idea ? 
>>>> 
>>>> I would consider using GMG. As Mark says, AMG is very fiddly with stretched grids. For Poisson, GMG works great and you seem to have regular grids.
>>>> 
>>>>   Thanks,
>>>> 
>>>>     Matt 
>>>> Thanks for your help, 
>>>> 
>>>> 
>>>> 
>>>> Victoria
>>>> 
>>>> 
>>>> 
>>>> On 7/10/20 12:48, Matthew Knepley wrote:
>>>>> On Wed, Oct 7, 2020 at 6:40 AM Victoria Hamtiaux <victoria.hamtiaux at uclouvain.be <mailto:victoria.hamtiaux at uclouvain.be>> wrote:
>>>>> Dear all,
>>>>> 
>>>>> 
>>>>> After the discretization of a poisson equation with purely Neumann (or 
>>>>> periodic) boundary conditions, I get a matrix which is singular.
>>>>> 
>>>>> 
>>>>> The way I am handling this is by using a NullSpace with the following 
>>>>> code :
>>>>> 
>>>>> MatNullSpace nullspace;
>>>>> MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace);
>>>>> MatSetNullSpace(p_solverp->A, nullspace);
>>>>> MatSetTransposeNullSpace(p_solverp->A, nullspace);
>>>>> MatNullSpaceDestroy(&nullspace);
>>>>> 
>>>>> 
>>>>> Note that I am using the hypre preconditionner BOOMERANG and the default 
>>>>> solver GMRES.
>>>>> 
>>>>> 
>>>>>      KSPCreate(PETSC_COMM_WORLD,&p_solverp->ksp);
>>>>>      KSPSetOperators(p_solverp->ksp, p_solverp->A, p_solverp->A);
>>>>>      PC prec;
>>>>>      KSPGetPC(p_solverp->ksp, &prec);
>>>>>      PCSetType(prec,PCHYPRE);//PCHYPRE seems the best
>>>>>      PCHYPRESetType(prec,"boomeramg"); //boomeramg is the best
>>>>>      KSPSetInitialGuessNonzero(p_solverp->ksp,PETSC_TRUE);
>>>>>      KSPSetFromOptions(p_solverp->ksp);
>>>>>      KSPSetTolerances(p_solverp->ksp, 1.e-10, 1e-10, PETSC_DEFAULT, 
>>>>> PETSC_DEFAULT);
>>>>>      KSPSetReusePreconditioner(p_solverp->ksp,PETSC_TRUE);
>>>>>      KSPSetUseFischerGuess(p_solverp->ksp,1,5);
>>>>>      KSPGMRESSetPreAllocateVectors(p_solverp->ksp);
>>>>>      KSPSetUp(p_solverp->ksp);
>>>>> 
>>>>> 
>>>>> 
>>>>> And this works fine when my grid is uniform, so that my matrix is 
>>>>> symmetric.
>>>>> 
>>>>> 
>>>>> But when I stretch the grid near the boundary (my grid is then 
>>>>> non-uniform), it doesn't work properly anymore. I suspect that by 
>>>>> stretching the grid, my matrix is not symmetric anymore and that it 
>>>>> might cause a problem.
>>>>> 
>>>>> Symmetry is a property of the operator, so you should be symmetric on your
>>>>> stretched grid. If not, I think you have the discretization wrong. You can check
>>>>> symmetry using
>>>>> 
>>>>>   https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatIsSymmetric.html <https://eur03.safelinks.protection.outlook.com/?url=https%3A%2F%2Fwww.mcs.anl.gov%2Fpetsc%2Fpetsc-current%2Fdocs%2Fmanualpages%2FMat%2FMatIsSymmetric.html&data=02%7C01%7Cvictoria.hamtiaux%40uclouvain.be%7C11ef103559a445e09bea08d86ada0450%7C7ab090d4fa2e4ecfbc7c4127b4d582ec%7C0%7C0%7C637376831896868460&sdata=fTsO10YkMdghL%2BW%2FAxnZieUN5mfkwPgcfQJXe7q3Is8%3D&reserved=0>
>>>>> 
>>>>> Also, if you suspect your discretization, you should probably do an MMS test to
>>>>> verify that you discretization converges at the correct rate.
>>>>> 
>>>>>   Thanks,
>>>>> 
>>>>>      Matt
>>>>>  
>>>>> I tried fixing the solution at an arbitrary point, but sometimes doing 
>>>>> this, I get errors near that fixed point. I 've seen on the petsc-users 
>>>>> forum that you usually don't recommend to fix a point, but I don't 
>>>>> really know how to proceed differently.
>>>>> 
>>>>> 
>>>>> What would you recommend to solve this problem?
>>>>> 
>>>>> 
>>>>> Thanks for your help,
>>>>> 
>>>>> 
>>>>> Best regards,
>>>>> 
>>>>> 
>>>>> Victoria
>>>>> 
>>>>> 
>>>>> 
>>>>> 
>>>>> -- 
>>>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>>>> -- Norbert Wiener
>>>>> 
>>>>> https://www.cse.buffalo.edu/~knepley/ <https://eur03.safelinks.protection.outlook.com/?url=http:%2F%2Fwww.cse.buffalo.edu%2F~knepley%2F&data=02%7C01%7Cvictoria.hamtiaux%40uclouvain.be%7C11ef103559a445e09bea08d86ada0450%7C7ab090d4fa2e4ecfbc7c4127b4d582ec%7C0%7C0%7C637376831896868460&sdata=VcIYJ8I3ObuuntK6fQz8XYL%2BjrzmFcw01TzLjAcvQnA%3D&reserved=0>
>>>> 
>>>> 
>>>> -- 
>>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>>> -- Norbert Wiener
>>>> 
>>>> https://www.cse.buffalo.edu/~knepley/ <https://eur03.safelinks.protection.outlook.com/?url=http:%2F%2Fwww.cse.buffalo.edu%2F~knepley%2F&data=02%7C01%7Cvictoria.hamtiaux%40uclouvain.be%7C11ef103559a445e09bea08d86ada0450%7C7ab090d4fa2e4ecfbc7c4127b4d582ec%7C0%7C0%7C637376831896878459&sdata=LZGtk6DGsCTXuSa0AGOAxMbiU%2FT1dTHMzomo6NsaDEw%3D&reserved=0>
>> 

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20201008/8bca9039/attachment-0001.html>


More information about the petsc-users mailing list