[petsc-users] Solving a coupled multi-physics problem.

Smith, Barry F. bsmith at mcs.anl.gov
Tue Jul 2 10:31:40 CDT 2019


    Rahul,

      This should be fine. But be careful since the original DMDA has run for 3 degrees of freedom per grid point while the u-displacement requires only one you will need to create a new DMDA identical to the first but with 1 dof (don't worry they and the parallel vectors they create will have the same parallel layout). 

      What you propose is a splitting scheme which introduces additional errors due to the split off of the u-displacement. You could possibly also consider a fully implicit approach where the time-stepper/SNES see all 4 variables directly with no splitting. This would require additional (approximate) Jacobian terms that represent the coupling between the pressure/flow and u-displacement.  If you do this you could use the higher implicit DAE/ODE integrators in PETSc instead of writing your own.

   Barry


> On Jul 2, 2019, at 6:30 AM, Rahul Samala via petsc-users <petsc-users at mcs.anl.gov> wrote:
> 
> Hello PetSc users,
> 
> I request suggestions on how to do a coupled problem, in my case coupling between multiphase flow in porous media and deformation of porous media.
> I have successfully written and validated the flow code with 3 degrees of freedom (S-saturation, P-pressure, T-Temperature)
> which uses SNES solver inside a time loop.  I have to use KSP solve for the deformation problem which has 1 degree of freedom (u-displacement).
> For the coupled problem, the flow part now requires the displacement vector and the deformation part requires
> the pressure field.
> 
> I have written a pseudo-code below. 
> My idea is to use the second solution vector (displacement) in FormFunctionLocal belonging to SNES.
> Since the DM used by SNES can be used for only one problem, a clone is made to be used for KSP.
> Is this approach the correct way to solve this coupled problem.
> 
> typedef struct {    
>     PetscScalar s; // Saturation
>     PetscScalar p; // Pressure
>     PetscScalar T; // Temperature
> } Field;
> 
> typedef struct {    
>     PetscScalar u; // Displacement
> } Field2;
> 
> FormFunctionLocal(DMDALocalInfo *info, Field *sol, Field2 *sol2, Field *f, AppCtx *user)
> {
> ........................
> }
> 
> main()
> {
>     DM   da, newda; 
>     Vec  sol;
>     Vec  sol2;
>     SNES  snes;
>     KSP   ksp;
> 
>     SNESCreate(......, &snes);
>     DMDACreate1d(........,&da);
>     .
>     .
>     SNESSetDM(snes, da);
>     DMClone(da, &newda);
>     KSPCreate(.........,&ksp);
>     KSPSetDM(ksp, newda);
>     .
>     .
>     DMDASetUniformCoordinates(da,...........);
>     DMDASNESSetFunctionLocal(da,...........,FormFunctionLocal);
>     .
>     DMCreateGlobalVector(da, &sol);
>     DMCreateGlobalVector(newda, &sol2);    // This sol2 of type Field2 is passed to FormFunctionLocal.
>     .
>     do { // Time loop
>     
>     SNESSolve(snes,....);
>     .
>     .
>     KSPSolve(ksp,......);    
>     .
>     .
>     }while(criteria);  
> }
> 
> Thank you,
> Rahul.
> 
> 
> 
> 



More information about the petsc-users mailing list