[petsc-users] query DMPlexCreateSection
Matthew Knepley
knepley at gmail.com
Fri Jan 3 10:31:20 CST 2014
On Fri, Jan 3, 2014 at 1:37 AM, Dharmendar Reddy <dharmareddy84 at gmail.com>wrote:
> Hello,
> I did not fully understand the flow. Need more help.
> If i look at the flow in DMPlexCreateSectionInitial
> (
> http://www.mcs.anl.gov/petsc/petsc-current/src/dm/impls/plex/plex.c.html#DMPlexCreateSection
> )
>
> 6023: for (d = 0; d <= dim; ++d) {
> 6024: numDofTot[d] = 0;
> 6025: for (f = 0; f < numFields; ++f) numDofTot[d] +=
> numDof[f*(dim+1)+d];
> 6026: }
> 6027: PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);
> 6028: if (numFields > 0) {
> 6029: PetscSectionSetNumFields(*section, numFields);
> 6030: if (numComp) {
> 6031: for (f = 0; f < numFields; ++f) {
> 6032: PetscSectionSetFieldComponents(*section, f, numComp[f]);
> 6033: }
> 6034: }
> 6035: }
> 6036: DMPlexGetChart(dm, &pStart, &pEnd);
> 6037: PetscSectionSetChart(*section, pStart, pEnd);
> 6038: for (d = 0; d <= dim; ++d) {
> 6039: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
> 6040: for (p = pStart; p < pEnd; ++p) {
> 6041: for (f = 0; f < numFields; ++f) {
> 6042: PetscSectionSetFieldDof(*section, p, f, numDof[f*(dim+1)+d]);
> 6043: }
> 6044: PetscSectionSetDof(*section, p, numDofTot[d]);
> 6045: }
> 6046: }
>
> In line 6027, number of fields for entire section is defined. Now,
> dose the flow in your earlier replies requires me to do that ? I.e,
> define numFields to 2. Then
>
> Inside the loop 6040, i need to check if the point belongs to region 1 or
> 2.
>
> If it belongs to region 1 then numFields is 1 and if it belongs to
> region 2, numFields is 2.
>
> Where do i set the number of fields in the flow you mentioned in the
> earlier email ?
>
The number of fields is always 2 and the chart is always the whole domain.
Its just that
one field only takes values on part of the chart. The other points will
have dof = 0.
Matt
> Thanks
> Reddy
>
> On Thu, Jan 2, 2014 at 9:56 AM, Matthew Knepley <knepley at gmail.com> wrote:
> > On Thu, Jan 2, 2014 at 4:11 AM, Dharmendar Reddy <
> dharmareddy84 at gmail.com>
> > wrote:
> >>
> >> Hello,
> >> I am trying to use DMPlexCreateSection from fortran. I was
> >> able to solve a Poisson equation earlier using DMPlex for handling
> >> mesh data.
> >>
> >> Now i need to solve a system of equations: (poisson + continuity
> >> equations)
> >>
> >> - div (grad phi) = (C + n) --- (1)
> >> div (J ) = 0 --------------- (2)
> >> J = n grad(phi)
> >> ( C is constant)
> >> Simulation domain is defined as, rectangular regions 1 and 2 shown
> below.
> >>
> >> -------------
> >> | 1 |
> >> -------------
> >> | |
> >> | 2 |
> >> | |
> >> --------------
> >>
> >> Now, equation 1 is defined in region 1 and 2
> >> and equation 2 is defined only for region 2.
> >>
> >> How do i setup the section ? DMPlexCreateSection applies the given
> >> DOF layout per cell to all cells in the mesh.
> >>
> >> I am not sure if all the function calls inside
> >> DMPlexCreateSectionIntial and DMPlexCreateSectionBCDof are accessible
> >> from Fotran.
> >
> >
> > You don't want them anyway since they apply to the whole domain. The
> control
> > flow could be:
> >
> > DMPlexSetChart()
> > <loop over region 1>
> > DMPlexSetDof() and DMPlexSetFieldDof()
> > <loop over region 2>
> > DMPlexSetDof() and DMPlexSetFieldDof()
> > <loop over BC>
> > DMPlexSetConstraintDof() and DMPlexSetFieldConstraintDof()
> > DMPlexSetUp()
> > <loop over BC>
> > DMPlexSetConstraintIndices() and DMPlexSetFieldConstraintIndices()
> >
> > We could try and package some of this up if it looks generic.
> >
> > Thanks,
> >
> > Matt
> >
> >>
> >> Thanks
> >> Reddy
> >>
> >> --
> >> -----------------------------------------------------
> >> Dharmendar Reddy Palle
> >
> >
> >
> >
> > --
> > 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
>
>
>
> --
> -----------------------------------------------------
>
--
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/20140103/7ea90c1f/attachment-0001.html>
More information about the petsc-users
mailing list