[petsc-users] PetscFEIntegrateBdResidual_Basic

Matthew Knepley knepley at gmail.com
Tue Nov 21 05:32:23 CST 2017


On Mon, Nov 20, 2017 at 7:19 PM, David Fuentes <fuentesdt at gmail.com> wrote:

> Sorry, If I call this AddBoundary Twice with different markers/ids, the
> bcFunc seems to be applied to set values for DM_BC_ESSENTIAL not
> DM_BC_NATURAL?
>
> Can I call    ierr = PetscDSSetBdResidual(prob, 0, f0_bd_u, f1_bd_zero)
> twice to set the f0 functions different for each boundary ?
>

Okay, this is the problem I was citing before. AddBoundary() does take
natural conditions, but it just ignores the inhomogeneous ones, which are
properly
integrated in the weak form. However, DS only takes a single pointfunction
and integrates over the whole boundary. Thus your fix. I will get this
fixed, but
since it entails some programming, it might take me a while. Sorry about
that.

  Thanks,

     Matt


>
> /*@C
>   PetscDSAddBoundary - Add a boundary condition to the model
>
>   Input Parameters:
> + ds          - The PetscDS object
> . type        - The type of condition, e.g. DM_BC_ESSENTIAL/DM_BC_ESSENTIAL_FIELD
> (Dirichlet), or DM_BC_NATURAL (Neumann)
> . name        - The BC name
> . labelname   - The label defining constrained points
> . field       - The field to constrain
> . numcomps    - The number of constrained field components
> . comps       - An array of constrained component numbers
> . bcFunc      - A pointwise function giving boundary values
> . numids      - The number of DMLabel ids for constrained points
> . ids         - An array of ids for constrained points
> - ctx         - An optional user context for bcFunc
>
>   Options Database Keys:
> + -bc_<boundary name> <num> - Overrides the boundary ids
> - -bc_<boundary name>_comp <num> - Overrides the boundary components
>
>   Level: developer
>
> .seealso: PetscDSGetBoundary()
> @*/
>
> On Mon, Nov 20, 2017 at 11:33 AM, David Fuentes <fuentesdt at gmail.com>
> wrote:
>
>> Thanks! will give that a try!
>> df
>>
>> On Mon, Nov 20, 2017 at 11:23 AM, Matthew Knepley <knepley at gmail.com>
>> wrote:
>>
>>> On Mon, Nov 20, 2017 at 12:20 PM, David Fuentes <fuentesdt at gmail.com>
>>> wrote:
>>>
>>>> Thanks for the quick reply. Indeed it does work like this. I have added
>>>> a location dependence on the boundary to differentiate the two. However,
>>>> when my mesh moves then the BC will be applied incorrectly.
>>>>
>>>
>>> I don't think you should have to do that. You can call AddBoundary
>>> twice. Once with a function with the first branch and marker value 3, and
>>> the next with the second branch and marker value 2. Does
>>> that make sense?
>>>
>>>   Thanks,
>>>
>>>      Matt
>>>
>>>> *//  PetscFEIntegrateBdResidual_Basic DMPlexComputeBdResidual_Internal*
>>>>
>>>>
>>>> *static* *void** f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,*
>>>>
>>>>                     *const** PetscInt uOff[], **const** PetscInt
>>>> uOff_x[], **const** PetscScalar u[], **const** PetscScalar u_t[], *
>>>> *const** PetscScalar u_x[],*
>>>>
>>>>                     *const** PetscInt aOff[], **const** PetscInt
>>>> aOff_x[], **const** PetscScalar a[], **const** PetscScalar a_t[], *
>>>> *const** PetscScalar a_x[],*
>>>>
>>>> *                    PetscReal t, **const** PetscReal x[], **const**
>>>> PetscReal n[], PetscInt numConstants, **const** PetscScalar
>>>> constants[], PetscScalar f0[])*
>>>>
>>>> *{ *
>>>>
>>>> *  PetscInt d;*
>>>>
>>>>   *double** radius = **0.0**;*
>>>>
>>>>   *const* *double**  zthresh = (**3.** - **4.5**/**2.**)***.01**; **//
>>>> [m]*
>>>>
>>>>   *for** (d = **0**; d < **2**; ++d) radius += x[d]*x[d];*
>>>>
>>>> *  radius = sqrt(radius); *
>>>>
>>>>   *if** ( radius > **.001**/**2.**  ) *
>>>>
>>>> *    {f0[**0**] =  constants[**2**] * u[**0**];}*
>>>>
>>>>   *else*
>>>>
>>>> *    {*
>>>>
>>>> *     f0[**0**] = -constants[**3**]       ; *
>>>>
>>>> *    }*
>>>>
>>>> *} *
>>>>
>>>>
>>>> On Mon, Nov 20, 2017 at 11:02 AM, Matthew Knepley <knepley at gmail.com>
>>>> wrote:
>>>>
>>>>> On Mon, Nov 20, 2017 at 11:46 AM, David Fuentes <fuentesdt at gmail.com>
>>>>> wrote:
>>>>>
>>>>>> *Is there a way to pass the boundary set id to the function pointers
>>>>>> for the residual evaluation on the boundary ?*
>>>>>>
>>>>>> *https://bitbucket.org/petsc/petsc/src/d89bd21cf2b5366df29efb6006298d2bc22fb509/src/dm/dt/interface/dtfe.c?at=master&fileviewer=file-view-default#dtfe.c-4245
>>>>>> <https://bitbucket.org/petsc/petsc/src/d89bd21cf2b5366df29efb6006298d2bc22fb509/src/dm/dt/interface/dtfe.c?at=master&fileviewer=file-view-default#dtfe.c-4245>*
>>>>>>
>>>>>> *I want to pass the boundary condition/constraint ID (ids):
>>>>>> PetscErrorCode PetscDSAddBoundary(PetscDS ds, DMBoundaryConditionType type,
>>>>>> const char name[], const char labelname[], PetscInt field, PetscInt
>>>>>> numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids,
>>>>>> const PetscInt *ids, void *ctx)*
>>>>>>
>>>>>> *to the functions for the residual evaluation on the boundary.*
>>>>>>
>>>>>>
>>>>>> *For example, I have two side sets in an exodus file. I want to
>>>>>> implement Neumann boundary conditions on side set = 2 and Mixed/Cauchy BC
>>>>>> on side set = 3. Or similarly use different*
>>>>>>
>>>>>> *gmsh BC tags for Neumann/Mixed BC.*
>>>>>>
>>>>> I am not completely against this, but let me respond with my rationale
>>>>> first. What I thought you would do, is call AddBoundary() twice. Once with
>>>>> the
>>>>> Neumann function and value 2, and once with the Cauchy function and
>>>>> value 3. Does that not work in your situation?
>>>>>
>>>>> Also, I am refectoring this right now because a DS object can only
>>>>> take a single boundary integral point function (which is a pain for
>>>>> inhomogeneous Neumann),
>>>>> so I welcome input.
>>>>>
>>>>>   Thanks,
>>>>>
>>>>>     Matt
>>>>>
>>>>> --
>>>>> 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
>>>>>
>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>> <http://www.caam.rice.edu/~mk51/>
>>>>>
>>>>
>>>>
>>>
>>>
>>> --
>>> 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
>>>
>>> https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
>>>
>>
>>
>


-- 
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

https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171121/3d52e148/attachment.html>


More information about the petsc-users mailing list