[petsc-users] Using DMPlex for Darcys law

Matthew Knepley knepley at gmail.com
Wed Aug 6 19:27:22 CDT 2014


On Wed, Aug 6, 2014 at 5:55 PM, Justin Chang <jychang48 at gmail.com> wrote:

> Hi,
>
> So in ex62 of SNES it solves the stokes equation with P2-P1 elements. I
> want to solve the darcy equation with P2-P1 elements in a similar fashion.
>
> The difference between stokes and darcy is that in the weak form of darcy,
> the velocity-velocity term is v*u whereas in stokes it's grad[v]*grad[u]
>
> In ex62, the function g3_uu (lines 172-183) formulates the grad[v]*grad[u]
> term but how would I do it for Darcy's equation where its simply v*u
> (assuming drag coefficient is scalar quantity of 1.0)?
>

So you would want to define g0 instead of g3, since its the u-u term:

void g0_uu(const PetscScalar u[], const PetscScalar u_t[], const
PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const
PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])
{
  const PetscInt Ncomp = spatialDim;
  PetscInt       compI;

  for (compI = 0; compI < Ncomp; ++compI) g0[compI] = u[compI];
}

and then register it in SetupProblem:

  ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL,  NULL,
 NULL);CHKERRQ(ierr);

and do not forget to change the residual definition in f0/f1.

  Thanks,

     Matt


> Thanks,
> Justin
>



-- 
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
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140806/74689d1e/attachment.html>


More information about the petsc-users mailing list