[petsc-users] Multigrid Preconditioned CG on cell-centered Poisson with DMStag
Barry Smith
bsmith at petsc.dev
Wed Oct 2 09:42:11 CDT 2024
Using
PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dmPressure), PETSC_TRUE, 0, nullptr, &matNullSpace));
should be fine (the other approach is equivalent).
Any change if you use -ksp_type gmres?
Please send all the output from a sample run with -ksp_monitor_true_residual -ksp_view
> On Oct 1, 2024, at 6:42 PM, Gautam Luhana <gkluhana at cs.ubc.ca> wrote:
>
> I did think that I accounted for that but now I see that I am getting the same number of iterations in the pressure space with and without the call to MatSetNullspace so maybe I am not doing that correctly and that's where the problem is.
>
> I tried it in the following way first:
>
> PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dmPressure), PETSC_TRUE, 0, nullptr, &matNullSpace));
>
> PetscCall(MatSetNullSpace(*A, matNullSpace));
>
> and then also as:
>
> PetscCall(DMGetGlobalVector(dmPressure, &constantPressure));
> PetscCall(VecSet(constantPressure, 1.0));
> PetscCall(VecNorm(constantPressure, NORM_2, &nrm));
> PetscCall(VecScale(constantPressure, 1.0 / nrm));
> PetscCall(DMCreateGlobalVector(dmPressure, &basis));
> PetscCall(VecCopy(constantPressure, basis));
> PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dmPressure), PETSC_FALSE, 1, &basis, &matNullSpace));
> PetscCall(MatSetNullSpace(*A, matNullSpace));
>
> Is this the right way? are there other options that need to be set?
>
> -Gautam
>
> On 2024-10-01 12:47 p.m., Barry Smith wrote:
>> [CAUTION: Non-UBC Email]
>>
>> Are you handling the null space? Like with MatSetNullSpace()?
>>
>> Can you try src/ksp/ksp/tutorials/ex32.c it solves the same problem but with a simpler DMDA. I don't see any reason why using DMStag with the cell-centered values should be different.
>>
>> Barry
>>
>>
>>> On Oct 1, 2024, at 3:37 PM, Gautam Luhana <gkluhana at cs.ubc.ca> wrote:
>>>
>>> Hello,
>>>
>>> I am solving the stokes equations on the MAC grid using DMStag with an approximate Schur complement based preconditioner. The preconditioner needs a subdomain solve of the Poisson equation on both the velocity as well as the pressure space. I am using multigrid preconditioned CG for both the subdomain solves. While the subdomain solve on the velocity space converges in 4-5 iterations and is independent of mesh-size, the pressure solve is not mesh-size independent. Relevant information:
>>>
>>> Pressure space operator has Neumann boundary conditions.
>>>
>>> Coarse-grid operators are being generated using the Galerkin approach (A = P*Af*Pt) for both of the subdomains.
>>>
>>> For the pressure space, the jacobi and sor smoothers were tried with several omega and damping parameter choices.
>>>
>>> Symmetric smoothing (same number of pre- and post-smoothing iterations)
>>>
>>> Has anyone tried solving the poisson operator on the cell-centered pressure space with DMStag? Or know what might be going wrong with the mg-cg convergence?
>>>
>>> Thanks!
>>>
>>> -Gautam
>>>
>>>
More information about the petsc-users
mailing list