Jed, Junchao, Thanks for getting me back!
Junchao, I cannot provide you with an input file as I am plugging just part of that code into mine. I am building blocks of the full matrix on my own, not reading from the file unfortunately. I can share a branch of my code with a 3x3 cavity toy case, but that would require compiling it and it might be a bit inconvenient for you. Jed, I am using fieldsplit indeed, it's just the matrix that is assembled in a 2x2 block form. Here is the full backtrace: *[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------* *[0]PETSC ERROR: No support for this operation for this object type* *[0]PETSC ERROR: No method getinfo for Mat of type nest* *[0]PETSC ERROR: WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_p_pc_type value: hypre source: command line* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_u_pc_type value: bjacobi source: command line* *[0]PETSC ERROR: See https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ <https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ > for trouble shooting.* *[0]PETSC ERROR: Petsc Release Version 3.20.4, unknown * *[0]PETSC ERROR: flubio_coupled on a linux_x86_64 named localhost.localdomain by edo Thu Oct 3 14:17:20 2024* *[0]PETSC ERROR: Configure options PETSC_ARCH=linux_x86_64 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_repo/openmpi-5.0.1/build/* *[0]PETSC ERROR: #1 MatGetInfo() at /home/edo/software_repo/petsc/src/mat/interface/matrix.c:3005* *[0]PETSC ERROR: #2 PCSetUp_GAMG() at /home/edo/software_repo/petsc/src/ksp/pc/impls/gamg/gamg.c:619* *[0]PETSC ERROR: #3 PCSetUp() at /home/edo/software_repo/petsc/src/ksp/pc/interface/precon.c:1080* *[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------* *[0]PETSC ERROR: It appears a new error in the code was triggered after a previous error, possibly because:* *[0]PETSC ERROR: - The first error was not properly handled via (for example) the use of* *[0]PETSC ERROR: PetscCall(TheFunctionThatErrors()); or* *[0]PETSC ERROR: - The second error was triggered while handling the first error.* *[0]PETSC ERROR: Above is the traceback for the previous unhandled error, below the traceback for the next error* *[0]PETSC ERROR: ALL ERRORS in the PETSc libraries are fatal, you should add the appropriate error checking to the code* *[0]PETSC ERROR: No support for this operation for this object type* *[0]PETSC ERROR: No method getinfo for Mat of type nest* *[0]PETSC ERROR: WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_p_pc_type value: hypre source: command line* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_u_pc_type value: bjacobi source: command line* *[0]PETSC ERROR: See https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ <https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ > for trouble shooting.* *[0]PETSC ERROR: Petsc Release Version 3.20.4, unknown * *[0]PETSC ERROR: flubio_coupled on a linux_x86_64 named localhost.localdomain by edo Thu Oct 3 14:17:20 2024* *[0]PETSC ERROR: Configure options PETSC_ARCH=linux_x86_64 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_repo/openmpi-5.0.1/build/* *[0]PETSC ERROR: #1 MatGetInfo() at /home/edo/software_repo/petsc/src/mat/interface/matrix.c:3005* *[0]PETSC ERROR: #2 PCSetUp_GAMG() at /home/edo/software_repo/petsc/src/ksp/pc/impls/gamg/gamg.c:619* *[0]PETSC ERROR: #3 PCSetUp() at /home/edo/software_repo/petsc/src/ksp/pc/interface/precon.c:1080* *[0]PETSC ERROR: #4 KSPSetUp() at /home/edo/software_repo/petsc/src/ksp/ksp/interface/itfunc.c:415* *[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------* *[0]PETSC ERROR: It appears a new error in the code was triggered after a previous error, possibly because:* *[0]PETSC ERROR: - The first error was not properly handled via (for example) the use of* *[0]PETSC ERROR: PetscCall(TheFunctionThatErrors()); or* *[0]PETSC ERROR: - The second error was triggered while handling the first error.* *[0]PETSC ERROR: Above is the traceback for the previous unhandled error, below the traceback for the next error* *[0]PETSC ERROR: ALL ERRORS in the PETSc libraries are fatal, you should add the appropriate error checking to the code* *[0]PETSC ERROR: No support for this operation for this object type* *[0]PETSC ERROR: No method getinfo for Mat of type nest* *[0]PETSC ERROR: WARNING! There are unused option(s) set! Could be the program crashed before usage or a spelling mistake, etc!* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_p_pc_type value: hypre source: command line* *[0]PETSC ERROR: Option left: name:-UPeqn_fieldsplit_u_pc_type value: bjacobi source: command line* *[0]PETSC ERROR: See https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ <https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ > for trouble shooting.* *[0]PETSC ERROR: Petsc Release Version 3.20.4, unknown * *[0]PETSC ERROR: flubio_coupled on a linux_x86_64 named localhost.localdomain by edo Thu Oct 3 14:17:20 2024* *[0]PETSC ERROR: Configure options PETSC_ARCH=linux_x86_64 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_repo/openmpi-5.0.1/build/* *[0]PETSC ERROR: #1 MatGetInfo() at /home/edo/software_repo/petsc/src/mat/interface/matrix.c:3005* *[0]PETSC ERROR: #2 PCSetUp_GAMG() at /home/edo/software_repo/petsc/src/ksp/pc/impls/gamg/gamg.c:619* *[0]PETSC ERROR: #3 PCSetUp() at /home/edo/software_repo/petsc/src/ksp/pc/interface/precon.c:1080* *[0]PETSC ERROR: #4 KSPSetUp() at /home/edo/software_repo/petsc/src/ksp/ksp/interface/itfunc.c:415* *[0]PETSC ERROR: #5 KSPSolve_Private() at /home/edo/software_repo/petsc/src/ksp/ksp/interface/itfunc.c:832* *[0]PETSC ERROR: #6 KSPSolve() at /home/edo/software_repo/petsc/src/ksp/ksp/interface/itfunc.c:1079* This the code I am doing right now: * M = this%setMomentumBlockMatrix() C = this%setPressureBlockMatrix() G = this%setPressureGradBlockMatrix() D = this%setDivBlockMatrix() bhat = this%setPressureMomentumRhs() !------------------! ! Create Ahat ! !-----------------! ! Bottom left block of Ahat Mat_array(4) = M ! Extraction of the diagonal of M call MatCreateVecs(M, PETSC_NULL_VEC, v, ierr) !v has the parallel distribution of M call MatGetDiagonal(M, v, ierr) call MatCreateDiagonal(v, diag_2M, ierr) call MatScale(diag_2M, 2.0, ierr) ! store 2*diagonal part of M call MatConvert(diag_2M, MATAIJ, MAT_INPLACE_MATRIX, diag_2M, ierr) call VecReciprocal(v, ierr) ! Creation of D_M_inv_G = D_M_inv*G call MatDuplicate(G,MAT_COPY_VALUES, D_M_inv_G, ierr) ! D_M_inv_G contains G call MatCreateVecs(D_M_inv_G, PETSC_NULL_VEC, v_redistributed, ierr) ! v_redistributed has the parallel distribution of D_M_inv_G call VecGetOwnershipRange(v, col_min, col_max, ierr) call ISCreateStride(PETSC_COMM_WORLD, col_max-col_min, col_min, 1, is_from, ierr) call VecGetOwnershipRange(v_redistributed, col_min, col_max, ierr) call ISCreateStride(PETSC_COMM_WORLD, col_max-col_min, col_min, 1, is_to, ierr) call VecScatterCreate(v,is_from,v_redistributed,is_to, scat, ierr) call VecScatterBegin(scat, v, v_redistributed,INSERT_VALUES,SCATTER_FORWARD, ierr) call VecScatterEnd( scat, v, v_redistributed,INSERT_VALUES,SCATTER_FORWARD, ierr) call MatDiagonalScale( D_M_inv_G, v_redistributed, PETSC_NULL_VEC, ierr) ! D_M_inv_G contains D_M_inv*G call ISDestroy(is_to, ierr) call ISDestroy(is_from, ierr) ! Creation of Chat (PETSC_DEFAULT_REAL -> PETSC_DETERMINE) call MatMatMult(D, D_M_inv_G, MAT_INITIAL_MATRIX, PETSC_DEFAULT_REAL, Chat, ierr) ! Chat contains D*D_M_inv*G call MatAXPY( Chat, 1.0, C, SUBSET_NONZERO_PATTERN, ierr) ! Chat contains C + D*D_M_inv*G Mat_array(1) = Chat ! Top left block of Ahat ! Creation of Ghat call MatMatMult(M, D_M_inv_G, MAT_INITIAL_MATRIX, PETSC_DEFAULT_REAL, Ghat, ierr) ! Ghat contains M*D_M_inv*G call MatAYPX( Ghat, -1.0, G, UNKNOWN_NONZERO_PATTERN, ierr) ! Ghat contains G - M*D_M_inv*G Mat_array(3) = Ghat ! Bottom left block of Ahat ! Creation of -D call MatScale(D,-1.0, ierr) Mat_array(2) = D ! Top right block of Ahat ! Creation of Ahat = transformed+ reordered A_input call MatCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL_IS, 2, PETSC_NULL_IS, Mat_array, Ahat, ierr) ! Creation of Pmat Mat_array(4) = diag_2M Mat_array(2) = PETSC_NULL_MAT ! Cancel top right block call MatCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL_IS, 2, PETSC_NULL_IS, Mat_array, Pmat, ierr); ! Finalisation of the preconditioner call ISCreateStride(PETSC_COMM_WORLD, numberOfElements, 0, 1, is_P_hat, ierr); call ISCreateStride(PETSC_COMM_WORLD, (3-bdim)*numberOfElements, numberOfElements, 1, is_U_hat, ierr); !----------------! ! KSP PART ! !----------------! call KSPSetOperators(this%ksp, Ahat, Pmat, ierr) call KSPGetPC(this%ksp, pc, ierr) call PCSetFromOptions(pc, ierr) call PCSetUp(pc, ierr) call KSPSetFromOptions(this%ksp, ierr) call KSPSetUp(this%ksp, ierr) call VecDuplicate(bhat, xhat, ierr) ! I think it makes sense to move this part where fieldsplit preconditioner is defined call PCFieldSplitSetIS(pc, "phat", is_P_hat, ierr) call PCFieldSplitSetIS(pc, "uhat", is_U_hat, ierr) call KSPSolve(this%ksp, bhat, xhat, ierr)* KSP id FGMRES + fieldsplit with default settings. Many thanks! Il giorno gio 3 ott 2024 alle ore 19:01 Junchao Zhang < junchao.zh...@gmail.com> ha scritto: > Could you also provide input files to run your code? > > --Junchao Zhang > > > On Thu, Oct 3, 2024 at 5:14 AM Edoardo alinovi <edoardo.alin...@gmail.com> > wrote: > >> Hi Barry, >> >> I am trying to port this preconditioner I found presented at the last >> PETSc conference into may code: >> >> https://urldefense.us/v3/__https://github.com/ndjinga/SOLVERLAB/blob/master/CoreFlows/examples/C/SaddlePointPreconditioner/SaddlePointLinearSolver.c__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF2QaZ-jc$ >> >> <https://urldefense.us/v3/__https://github.com/ndjinga/SOLVERLAB/blob/master/CoreFlows/examples/C/SaddlePointPreconditioner/SaddlePointLinearSolver.c__;!!G_uCfscf7eWS!a-JnrBsta75n_WDtxUNd7Xb8_m7tg1pWSxwWNTmUBSfpBwISq4IxJkbM6QgADd-x23Q_2s-d6d6bOWCKGRMxMzHyc9DDm7A$> >> >> So far so good, I am just translating C to Fortran, however when I try to >> solve my linear system with: >> >> call KSPSetOperators(ksp, Ahat, Pmat, ierr) >> call KSPSolve(ksp, bhat, xhat, ierr) >> >> with Ahat and Pmat being indeed two nest MatNest, created from an array >> of Mat as: >> >> call MatCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL_IS, 2, PETSC_NULL_IS, >> Mat_array, Ahat, ierr) >> call MatCreateNest(PETSC_COMM_WORLD, 2, PETSC_NULL_IS, 2, PETSC_NULL_IS, >> Mat_array, Pmat, ierr) >> >> >> However I am getting this error on getinfo: >> >> from C to Fortran >> >> [0]PETSC ERROR: --------------------- Error Message >> -------------------------------------------------------------- >> [0]PETSC ERROR: No support for this operation for this object type >> [0]PETSC ERROR: No method getinfo for Mat of type nest >> [0]PETSC ERROR: See >> https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a6ejiT9ML4BpeAAcqVuIMX-J7sA75OiyYcHeUd6-INXK8sR-gcQPA3ndqsgpzphK4tIBnE13dbzycfmUkwST7aXF6oxXipY$ >> >> <https://urldefense.us/v3/__https://petsc.org/release/faq/__;!!G_uCfscf7eWS!a-JnrBsta75n_WDtxUNd7Xb8_m7tg1pWSxwWNTmUBSfpBwISq4IxJkbM6QgADd-x23Q_2s-d6d6bOWCKGRMxMzHyRUMnyvE$> >> for trouble shooting. >> [0]PETSC ERROR: Petsc Release Version 3.20.4, unknown >> >> How can I get rid of this issue? >> >> Many thanks! >> >>