Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Update to [email protected] #126

Merged
merged 1 commit into from
Jan 6, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@ Random = "9a3f8284-a2c9-5f02-9a11-845980a1fd5c"
SpecialFunctions = "276daf66-3868-5448-9aa4-cd146d93841b"

[compat]
Gridap = "0.16.4"
Gridap = "0.17"
GridapGmsh = "0.4"
SpecialFunctions = "1"
julia = "1.3"
Expand Down
36 changes: 15 additions & 21 deletions src/fsi_tutorial.jl
Original file line number Diff line number Diff line change
Expand Up @@ -176,9 +176,18 @@ const μ_f = 0.01
# Q \doteq \{ q \in C^0(\Omega):\ q|_T\in P_{k-1}(T) \text{ for all } T\in\mathcal{T}_{\rm F}\}.
# ```

# Before creating the FE spaces, we first need to define the discrete model of each of the sub-domains, which are constructed restricting the global discrete model to the corresponding part. This is done by calling the `DiscreteModel` function with the desired geometrical part label, i.e. `"solid"` and `"fluid"`, respectively.
model_solid = DiscreteModel(model,tags="solid")
model_fluid = DiscreteModel(model,tags="fluid")
# Before creating the FE spaces, we first need to define the triangulation of each of the sub-domains, which are constructed restricting the global triangulation to the corresponding part. This is done by calling the `Triangulation` (or equivalently `Interior`) function with the desired geometrical part label, i.e. `"solid"` and `"fluid"`, respectively. Here we create the triangulation of the global domain, $\mathcal{T}$, and the solid and fluid triangulations, $\mathcal{T}_{\rm F}$ and $\mathcal{T}_{\rm S}$.
Ω = Interior(model)
Ω_s = Interior(model,tags="solid")
Ω_f = Interior(model,tags="fluid")

# We also generate the triangulation and associated outer normal field for the outlet boundary, $\Gamma_{\rm F,N_{out}}$, which will be used to impose a Neumann condition.
Γ_out = BoundaryTriangulation(model,tags="outlet")
n_Γout = get_normal_vector(Γ_out)

# Finally, to impose the interface condition between solid and fluid, we will need the triangulation and normal field of such interface, $\Gamma_{\rm FS}$. The interface triangulation is generated by calling the `InterfaceTriangulation` function specifying the two interfacing domain triangulations. Note that the normal field will point outwards with respect to the first entry.
Γ_fs = InterfaceTriangulation(Ω_f,Ω_s)
n_Γfs = get_normal_vector(Γ_fs)

# In what follows we will assume a second-order veloticty interpolation, i.e. $k=2$
k = 2
Expand All @@ -189,20 +198,20 @@ reffeₚ = ReferenceFE(lagrangian,Float64,k-1)

# Having set up all the ingredients, we can create the different FE spaces for the test functions. For the velocity FE spaces we call the `TestFESpace` function with the corresponding discrete model, using the velocity reference FE `reffeᵤ` and conformity `:H1`. Note that we assign different Dirichlet boundary labels for the two different parts, generating the variational spaces with homogeneous Dirichlet boundary conditions, $V_{\rm F,0}$ and $V_{\rm S,0}$ .
Vf = TestFESpace(
model_fluid,
Ω_f,
reffeᵤ,
conformity=:H1,
dirichlet_tags=["inlet", "noSlip", "cylinder"])

Vs = TestFESpace(
model_solid,
Ω_s,
reffeᵤ,
conformity=:H1,
dirichlet_tags=["fixed"])

# For the pressure test FE space, we use the fluid discrete model, the pressure reference FE `reffeₚ` and `:C0` conformity.
Qf = TestFESpace(
model_fluid,
Ω_f,
reffeₚ,
conformity=:C0)

Expand All @@ -219,21 +228,6 @@ X = MultiFieldFESpace([Us,Uf,Pf])
# <a name="integration"></a>
# ```
# ### Numerical integration
# To define the quadrature rules used in the numerical integration of the different terms, we first need to generate the domain triangulation. Here we create the triangulation of the global domain, $\mathcal{T}$.
Ω = Triangulation(model)

# The solid and fluid triangulations, $\mathcal{T}_{\rm F}$ and $\mathcal{T}_{\rm S}$, are constructed from the discrete models restricted to the respective parts.
Ω_s = Triangulation(model_solid)
Ω_f = Triangulation(model_fluid)

# We also generate the triangulation and associated outer normal field for the outlet boundary, $\Gamma_{\rm F,N_{out}}$, which will be used to impose a Neumann condition.
Γ_out = BoundaryTriangulation(model,tags="outlet")
n_Γout = get_normal_vector(Γ_out)

# Finally, to impose the interface condition between solid and fluid, we will need the triangulation and normal field of such interface, $\Gamma_{\rm FS}$. The interface triangulation is generated by calling the `InterfaceTriangulation` function specifying the two interfacing domain models. Note that the normal field will point outwards with respect to the first entry.
Γ_fs = InterfaceTriangulation(model_fluid,model_solid)
n_Γfs = get_normal_vector(Γ_fs)

# Once we have all the triangulations, we can generate the quadrature rules to be applied each domain. This will be generated by calling the `Measure` function, that given a triangulation and an integration degree, it returns the Lebesgue integral measure $d\Omega$.
degree = 2*k
dΩₛ = Measure(Ω_s,degree)
Expand Down