[petsc-dev] DMComplex and Coordinates
Chris Eldred
chris.eldred at gmail.com
Wed Sep 26 16:47:23 CDT 2012
Looks great- thanks!
On Wed, Sep 26, 2012 at 3:46 PM, Matthew Knepley <knepley at gmail.com> wrote:
> On Wed, Sep 26, 2012 at 5:37 PM, Chris Eldred <chris.eldred at gmail.com>
> wrote:
>>
>> Will there be a way to set coordinates if I am using a homemade mesh
>> generator?
>
>
> Sure, here is me doing this by hand (from complexcreate.c):
>
> /* Build coordinates */
> ierr = DMComplexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
> ierr = PetscSectionSetChart(coordSection, numEdges, numEdges +
> numVertices);CHKERRQ(ierr);
> for (v = numEdges; v < numEdges+numVertices; ++v) {
> ierr = PetscSectionSetDof(coordSection, v, 2);CHKERRQ(ierr);
> }
> ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
> ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
> ierr = VecCreate(((PetscObject) dm)->comm, &coordinates);CHKERRQ(ierr);
> ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
> ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr);
> ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
> for (vy = 0; vy <= edges[1]; ++vy) {
> for (vx = 0; vx <= edges[0]; ++vx) {
> coords[(vy*(edges[0]+1)+vx)*2+0] = lower[0] + ((upper[0] -
> lower[0])/edges[0])*vx;
> coords[(vy*(edges[0]+1)+vx)*2+1] = lower[1] + ((upper[1] -
> lower[1])/edges[1])*vy;
> }
> }
> ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
> ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
>
>>
>> What about having coordinates for the edges/faces/cells as well (to be
>> interpreted as we want)?
>
>
> Yep, just give them a dof in the Section and set the values.
>
> Matt
>
>>
>> On Wed, Sep 26, 2012 at 3:34 PM, Matthew Knepley <knepley at gmail.com>
>> wrote:
>> > On Wed, Sep 26, 2012 at 5:24 PM, Chris Eldred <chris.eldred at gmail.com>
>> > wrote:
>> >>
>> >> On the subject:
>> >>
>> >> How do coordinates work for user-defined DMComplex objects? What
>> >> determines the Coordinate section (and which points does it contain?)?
>> >> How is the coordinate vector set?
>> >
>> >
>> > I am just rewriting this. Now it will be
>> >
>> > DMGetCoordinateDM(DM, DM *)
>> > DMGetCoordinates(DM, Vec *)
>> > DMGetCoordinatesLocal(DM, Vec *)
>> >
>> > and there is a convenience
>> >
>> > DMComplexGetCoordinateSection(DM, PetscSection *)
>> >
>> > By default, I read in the vertex coordinates from the mesh generator.
>> > However,
>> > the idea is that you can interpret the coordinates however you want. We
>> > just
>> > give you a DM and a vector. You can construct the section however you
>> > want,
>> > although without vertex coordinates some default things will not work,
>> > like
>> > VTK
>> > output.
>> >
>> > Matt
>> >
>> >>
>> >>
>> >> --
>> >> Chris Eldred
>> >> DOE Computational Science Graduate Fellow
>> >> Graduate Student, Atmospheric Science, Colorado State University
>> >> B.S. Applied Computational Physics, Carnegie Mellon University, 2009
>> >> chris.eldred at gmail.com
>> >
>> >
>> >
>> >
>> > --
>> > 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
>>
>>
>>
>> --
>> Chris Eldred
>> DOE Computational Science Graduate Fellow
>> Graduate Student, Atmospheric Science, Colorado State University
>> B.S. Applied Computational Physics, Carnegie Mellon University, 2009
>> chris.eldred at gmail.com
>
>
>
>
> --
> 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
--
Chris Eldred
DOE Computational Science Graduate Fellow
Graduate Student, Atmospheric Science, Colorado State University
B.S. Applied Computational Physics, Carnegie Mellon University, 2009
chris.eldred at gmail.com
More information about the petsc-dev
mailing list