[petsc-users] On PCFIELDSPLIT and its implementation

Matthew Knepley knepley at gmail.com
Wed Nov 9 07:12:22 CST 2022


On Wed, Nov 9, 2022 at 8:09 AM Edoardo alinovi <edoardo.alinovi at gmail.com>
wrote:

> Sure,
>
> I'll try on a 3x3 cavity. How can I print the ISs?
>

ISView() or PetscObjectViewFromOptions()

  Thanks,

     Matt


> Il Mer 9 Nov 2022, 14:07 Matthew Knepley <knepley at gmail.com> ha scritto:
>
>> On Wed, Nov 9, 2022 at 7:57 AM Edoardo alinovi <edoardo.alinovi at gmail.com>
>> wrote:
>>
>>> To be clear,
>>>
>>> You are suggesting to use ufields(0)=0, ufields(1)=1 and so on?
>>>
>>
>> I think you are right. Those should start from 1. However, your ISes do
>> not seem to cover
>> the whole matrix. Can you start with a very small problem so that you can
>> print them to
>> the screen?
>>
>>   Thanks,
>>
>>      Matt
>>
>>
>>> Il Mer 9 Nov 2022, 13:54 Edoardo alinovi <edoardo.alinovi at gmail.com> ha
>>> scritto:
>>>
>>>> Even in the fortran interface?
>>>>
>>>> Il Mer 9 Nov 2022, 13:52 Matthew Knepley <knepley at gmail.com> ha
>>>> scritto:
>>>>
>>>>> Fields are numbered from 0.
>>>>>
>>>>>   Thanks,
>>>>>
>>>>>      Matt
>>>>>
>>>>> On Wed, Nov 9, 2022 at 2:20 AM Edoardo alinovi <
>>>>> edoardo.alinovi at gmail.com> wrote:
>>>>>
>>>>>> Hello guys,
>>>>>>
>>>>>> I am getting this error while using fieldsplit:
>>>>>>
>>>>>> [3]PETSC ERROR: --------------------- Error Message
>>>>>> --------------------------------------------------------------
>>>>>>
>>>>>> *[3]PETSC ERROR: Nonconforming object sizes[3]PETSC ERROR: Local
>>>>>> column sizes 6132 do not add up to total number of columns 9200*
>>>>>> [3]PETSC ERROR: See https://petsc.org/release/faq/ for trouble
>>>>>> shooting.
>>>>>> [3]PETSC ERROR: Petsc Development GIT revision:
>>>>>> v3.18.1-191-g32ed6ae2ff2  GIT Date: 2022-11-08 12:22:17 -0500
>>>>>> [3]PETSC ERROR: flubio_coupled on a gnu named alienware by edo Wed
>>>>>> Nov  9 08:16:29 2022
>>>>>> [3]PETSC ERROR: Configure options PETSC_ARCH=gnu FOPTFLAGS=-O3
>>>>>> COPTFLAGS=-O3 CXXOPTFLAGS=-O3 -with-debugging=no -download-fblaslapack=1
>>>>>> -download-superlu_dist -download-mumps -download-hypre -download-metis
>>>>>> -download-parmetis -download-scalapack -download-ml -download-slepc
>>>>>> -download-hpddm -download-cmake
>>>>>> -with-mpi-dir=/home/edo/software/openmpi-4.1.1/build/
>>>>>> [3]PETSC ERROR: #1 MatCreateSubMatrix_MPIBAIJ_Private() at
>>>>>> /home/edo/software/petsc/src/mat/impls/baij/mpi/mpibaij.c:1987
>>>>>> [3]PETSC ERROR: #2 MatCreateSubMatrix_MPIBAIJ() at
>>>>>> /home/edo/software/petsc/src/mat/impls/baij/mpi/mpibaij.c:1911
>>>>>> [3]PETSC ERROR: #3 MatCreateSubMatrix() at
>>>>>> /home/edo/software/petsc/src/mat/interface/matrix.c:8340
>>>>>> [3]PETSC ERROR: #4 PCSetUp_FieldSplit() at
>>>>>> /home/edo/software/petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c:657
>>>>>> [3]PETSC ERROR: #5 PCSetUp() at
>>>>>> /home/edo/software/petsc/src/ksp/pc/interface/precon.c:994
>>>>>> [3]PETSC ERROR: #6 KSPSetUp() at
>>>>>> /home/edo/software/petsc/src/ksp/ksp/interface/itfunc.c:406
>>>>>> [3]PETSC ERROR: #7 KSPSolve_Private() at
>>>>>> /home/edo/software/petsc/src/ksp/ksp/interface/itfunc.c:825
>>>>>> [3]PETSC ERROR: #8 KSPSolve() at
>>>>>> /home/edo/software/petsc/src/ksp/ksp/interface/itfunc.c:1071
>>>>>>
>>>>>> Do you have any ideas? Probably something missing in my brief
>>>>>> implementation here:
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> *            call PCSetType(mypc, PCFIELDSPLIT, ierr)
>>>>>> call PCFieldSplitSetBlockSize(mypc, 4-bdim, ierr)           *
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> *            !2D, 3x3 block            if(bdim==1) then
>>>>>>   ufields(1) = 0                ufields(2) = 1                pfields(1) =
>>>>>> 2                call PCFieldSplitSetFields(mypc, "u", 2, ufields, ufields,
>>>>>> ierr)                call PCFieldSplitSetFields(mypc, "p", 1, pfields,
>>>>>> pfields, ierr)             ! 3D 4x4 block            else
>>>>>> ufields(1) = 0                ufields(2) = 1                ufields(3) = 2
>>>>>>               pfields(1) = 3                call
>>>>>> PCFieldSplitSetFields(mypc, "u", 3, ufields, ufields, ierr)
>>>>>> call PCFieldSplitSetFields(mypc, "p", 1, pfields, pfields, ierr)
>>>>>> endif                         ! Field split type ADDITIVE, MULTIPLICATIVE
>>>>>> (default), SYMMETRIC_MULTIPLICATIVE, SPECIAL, SCHUR            call
>>>>>> PCFieldSplitSetType(mypc, PC_COMPOSITE_SCHUR, ierr)*
>>>>>>
>>>>>> Thanks for the help!
>>>>>>
>>>>>
>>>>>
>>>>> --
>>>>> 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/>
>>>>>
>>>>
>>
>> --
>> 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/>
>>
>

-- 
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/20221109/fdf5eda9/attachment-0001.html>


More information about the petsc-users mailing list