Hello, My name is Roberto. I'm currently a professor in the Applied Mathematics Department of Univ. of São Paulo at São Carlos. I'm working on some fluid-structure interaction problems in which the solid will be modeled as a 1D beam. Our reference integration domain will be the interval [0,L]. We are planning to use firedrake and Hermite elements. The problem is nonlinear. I already have a working example. Values degrees of freedom are easy to impose by using DirichletBC(), however I'm not figuring out how to enforce derivatives degrees of freedom as boundary conditions (let say at x=0), which I need for a clamped beam. Any help or trick would be appreciated Regards -- Roberto F. Ausas Instituto de Ciências Matemáticas e de Computação - ICMC/USP São Carlos, SP - Brasil *http://www.lmacc.icmc.usp.br/~ausas* <http://www.lmacc.icmc.usp.br/~ausas>
Hi Roberto, I did this before. I projected some functions onto 1D Hermite spaces to see which dofs are the derivative value at the endpoints. In the end I implemented class FixHermiteDerivativeValue(DirichletBC): @utils.cached_property def nodes(self): if V.mesh().mpi_comm().size > 1: raise NotImplementedError return [1] # the derivative component of the left endpoint bcs = [FixHermiteDerivativeValue(V, 0, "on_boundary")] This only fixes the derivative at the left; I guess the derivative on the right will be [-1], but am not sure. The "on_boundary" is ignored; I just kept it to not have to override the constructor of DirichletBC. Hope this helps, Patrick On 06/04/2020 20:49, Roberto Federico Ausas wrote:
Hello,
My name is Roberto. I'm currently a professor in the Applied Mathematics Department of Univ. of São Paulo at São Carlos.
I'm working on some fluid-structure interaction problems in which the solid will be modeled as a 1D beam. Our reference integration domain will be the interval [0,L].
We are planning to use firedrake and Hermite elements. The problem is nonlinear.
I already have a working example. Values degrees of freedom are easy to impose by using DirichletBC(), however I'm not figuring out how to enforce derivatives degrees of freedom as boundary conditions (let say at x=0), which I need for a clamped beam.
Any help or trick would be appreciated Regards
-- Roberto F. Ausas Instituto de Ciências Matemáticas e de Computação- ICMC/USP São Carlos, SP - Brasil _http://www.lmacc.icmc.usp.br/~ausas_
_______________________________________________ firedrake mailing list firedrake@imperial.ac.uk https://mailman.ic.ac.uk/mailman/listinfo/firedrake
participants (2)
-
Patrick Farrell
-
Roberto Federico Ausas