[petsc-users] PetscFEIntegrateBdResidual_Basic
Matthew Knepley
knepley at gmail.com
Mon Nov 20 11:23:12 CST 2017
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/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171120/10a3e035/attachment-0001.html>
More information about the petsc-users
mailing list