[petsc-dev] How to use fieldsplit for a block system with sizes change dynamically in TS

Mebratu Wakeni Mebratu.Wakeni at glasgow.ac.uk
Tue Sep 22 11:20:25 CDT 2020


Hello,

I have the following issue regarding TS and field split.

In a TS setting, I have a linear problem with the IJacobian of the IFunction with a block form

| K   0 |
| 0   D |

Where K is sparse, while D is a diagonal. The sizes of the square matrices K and D are dynamic (varying between time steps). The fact that, K and D are completely decoupled, and D is just a diagonal motivated us to use fieldsplit techniques rather than a monolithic one. However, this needs to recreate a fieldsplit preconditioner at each time step using the indices of K. To achieve this, I use the TSSetPreStep(ts, fun) function to perform the fun(TS ts) routine which does following


1.     Get indices of K,



PetscInt *idx_k; // size_of_k



2.     Create an IS corresponding to K (IS is_k) using



 ISCreatGeneral(mpiComm, size_of_k, idx_k, PETSC_COPY_VALUES, &is_k);



3.     Then set a fieldsplit preconditioner PC pc using is_k through



PCFieldSplitSetIS(pc, NULL, is_k);



4.     Finally the is_k is properly destroyed through



ISDestroy(&is_lo);

With the above setup, I have also used the following runtime options

-ts_monitor
-ts_exact_final_time stepover
-ts_dt 0.1
-ts_max_time 1
-ts_adapt_type none

-ksp_type gmres
-pc_type fieldsplit
-pc_fieldsplit_diag_use_amat
-pc_fieldsplit_type additive
-fieldsplit_0_ksp_type preonly
-fieldsplit_0_pc_type lu
-fieldsplit_0_pc_factor_mat_solver_type mumps
-fieldsplit_1_ksp_type preonly
-fieldsplit_1_pc_type jacobi

-snes_atol 1e-8
-snes_rtol 1e-12
-ksp_monitor
-snes_monitor
-snes_lag_jacobian 1

Result: with the setup described above and the options, I expect to get ksp convergence in a single step (at every time step), however, what I am getting in actuality is that the fieldsplit works as expected only for the first timestep, after that it either take too many ksp steps to converge or diverge completely. I have a suspicion that, the way the fieldsplit set up at each time using the TSSetPreStep() somehow adds new split over previous splits. For this, I have also tried using TSSetPostStep() function to reset both the ksp and pc as follows

KSPRest(ksp);
PCReset(pc);

But these don’t help either.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20200922/a8f1b779/attachment.html>


More information about the petsc-dev mailing list