[petsc-users] Testing DMPlex new functionalities
Luc Berger-Vergiat
lb2653 at columbia.edu
Thu Jan 22 16:31:31 CST 2015
Ok, I think what I will do is the following:
void top_bc(const PetscReal x[], PetscScalar *u, void *ctx)
{
*u = 0.2;
}
void bottom_bc(const PetscReal x[], PetscScalar *u, void *ctx)
{
*u = 0.0;
}
ierr = DMPlexAddBoundary(dm, PETSC_TRUE,"bottom","Face Sets",0,(void
(*)())bottom_bc,1,ids[0],&user);CHKERRQ(ierr);
ierr = DMPlexAddBoundary(dm, PETSC_TRUE,"top","Face Sets",0,(void
(*)())top_bc,1,ids[1],&user);CHKERRQ(ierr);
Matt, have you started writing some documentation for these function?
If not I could probably get something started since I'm reading them.
Another question, how do I explain to my DMPlex that I am solving an
elasticity problem and hence my solution field will be a 2 components
vector and not a scalar? I haven't seen where you set this info in snes
example 12...
Best,
Luc
On 01/22/2015 05:19 PM, Matthew Knepley wrote:
> On Thu, Jan 22, 2015 at 4:17 PM, Luc Berger-Vergiat
> <lb2653 at columbia.edu <mailto:lb2653 at columbia.edu>> wrote:
>
> Yes you are right, I corrected for the machine epsilon.
> I also looked deeper in Petsc's source code and I think I
> understand better how you implemented the GMSH reader.
> For instance the label I really care about is "Face Sets" since it
> detects automatically the BC's for me.
> An issue is that "Face Sets" contains 4 strata and I only care
> about two of these strata.
> How can I get the names of the Strata and how can I loop over the
> strata to assign my BCs?
>
>
> Instead of "marker" give "Face Sets"
>
> Instead of 1, give 2
>
> Instead of &id give ids, where PetscInt ids[2] = {stratumValue1,
> stratumValue2}
>
> Thanks,
>
> Matt
>
> Best,
> Luc
>
> On 01/22/2015 05:12 PM, Matthew Knepley wrote:
>> On Thu, Jan 22, 2015 at 12:08 PM, Luc Berger-Vergiat
>> <lb2653 at columbia.edu <mailto:lb2653 at columbia.edu>> wrote:
>>
>> Hi all,
>> I'm trying to test some functionalities in DMPlex,
>> specifically I want to do a 2d fem simulation using a mesh
>> generated by GMSH.
>> At this point I am trying to understand how Dirichlet BC can
>> be applied.
>>
>> Can you explain to me what arguments I have to pass to:
>> DMPlexAddBoundary()?
>>
>>
>> ierr = DMPlexAddBoundary(dm, PETSC_TRUE,"wall","marker",(void
>> (*)())top_displacement,1,&id,user)
>>
>> PETSC_TRUE means that these unknowns will be eliminated from the
>> global system
>>
>> "wall" is the name of the BC
>>
>> "marker" is the name of the DMLabel that you want to use to
>> identify the mesh points
>>
>> 1 is the number of label IDs that identify constrained mesh points
>>
>> id is the label ID which identifies constrained mesh points
>>
>> top_displacement is the function used to set boundary values
>>
>> I am attaching a little bit of c code and my mesh.geo and
>> mesh.msh so that you understand the context.
>> The BCs I want to apply are displacement of 0.2 in the y
>> direction if y==1.0
>> 0.0 in the x direction if y==0.0
>>
>>
>> I might do this:
>>
>> void top_displacement(const PetscReal x[], PetscScalar *u, void *ctx)
>> {
>> *u = x[1] < PETSC_MACHINE_EPSILON ? 0.0 : 0.2;
>> }
>>
>> Thanks,
>>
>> Matt
>>
>>
>> --
>> Best,
>> Luc
>>
>>
>>
>>
>> --
>> 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/20150122/46ea3b00/attachment-0001.html>
More information about the petsc-users
mailing list