[petsc-users] Matrix assembly problem of overlapped DMPlex
Matthew Knepley
knepley at gmail.com
Sun Jul 23 15:28:10 CDT 2023
On Sun, Jul 23, 2023 at 11:54 AM 袁煕 <yuanxi at advancesoft.jp> wrote:
> 袁煕 <yuanxi at advancesoft.jp>
> 0:50 (0 分前)
> To PETSc
> Hi,
>
> I used PETSc to assemble a FEM stiff matrix of an overlapped (overlap=2)
> DMPlex and used the MUMPS solver to solve it. But I got a
> different solution by using 1 CPU and MPI parallel computation. I am
> wondering if I missed some necessary step or setting during my
> implementation.
>
> My calling process like follows
>
> 1. Generate matrix pkij from DMPlex dm_mesh
>
> ```
> PetscCallA( DMCreateMatrix(dm_mesh, pkij, ierr) )
> PetscCallA( MatSetBlockSize(pkij, 1, ierr))
> if( overlap>0 ) call MatSetOption(pkij,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE,ierr)
> PetscCallA( MatSetFromOptions(pkij, ierr))
>
> PetscCallA(MatGetLocalSize(pkij, m, n, ierr));
> PetscCallA(MatGetSize(pkij, gm, gn, ierr));
> PetscCallA(MatGetBlockSize(pkij, bs, ierr));
> ```
>
> 2. Do Matrix preallocation
>
> ```
> PetscCallA(MatCreate(PETSC_COMM_WORLD, preallocator, ierr));
> PetscCallA(MatSetType(preallocator, MATPREALLOCATOR, ierr));
> PetscCallA(MatSetSizes(preallocator, m, n, gm, gn, ierr));
> PetscCallA(MatSetBlockSize(preallocator, bs, ierr));
> PetscCallA(MatSetUp(preallocator, ierr));
>
> .......
>
> PetscCallA( MatPreallocatorPreallocate(preallocator, PETSC_TRUE, pkij, ierr) )
> PetscCallA( MatDestroy(preallocator,ierr) )
> ```
>
> 3. Assemble matrix pkij by calling MatSetValue of all overlapped elements.
>
> In my above implementation, I used MatSetOption(pkij,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE,ierr). Is that correct? Or even other options are needed?
>
>
Do you want to ignore off-process entries? Usually in FEM we share part of
the boundary with other processes, and thus ignoring the off-process
entries will not give the correct result.
Thanks,
Matt
>
>
> Thanks for your help.
>
> ------
>
> X. Yuan, Ph.D. of solid mechanics
>
> Advancesoft, Japan
>
> ---------
>
>
--
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/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230723/8c38e33c/attachment.html>
More information about the petsc-users
mailing list