I need them just on the faces that separate the materials (fluid and 
solid). The piece of code is the following, on an interface cell face I 
have:

cell->face(face_number)->get_dof_indices(local_face_dof_indices_fluid, 
fluid_fe_index); // select the fluid dofs sets from hp-dofs
cell->face(face_number)->get_dof_indices(local_face_dof_indices_solid, 
solid_fe_index); // select the solid dofs sets from hp-dofs 

// impose fluid_velocity = solid_velocity (where solid_velocity = 
(solid_displacement - solid_displacement_previous_time_step) / time_step)

for (unsigned int i = 0; i < local_face_dof_indices_fluid.size(); ++i)
{
if (fe_fluid->face_system_to_component_index(i).first < dim) // select just 
fluid velocity (the solution is stored as (fluid_velocity + fluid_pressure 
+ solid_displacement) with (dim + 1 + dim) components)
{
auto column_index = /*todo*/ ; // need to extract the right index to select 
the displacement value at the interface
constraints.add_line(local_face_dof_indices_fluid[i]);
constraints.add_entry(local_face_dof_indices_fluid[i], column_index, 1.0 / 
time_step);
constraints.set_inhomogeneity(local_face_dof_indices_fluid[i], - 
solution_previous_time_step(column_index) / time_step);
}
}

Basically I don't know how to extract the index (or indices) that point at 
the solid displacement value on the other side of the cell face. Is there a 
way to do that?

Thank you


Il giorno mercoledì 18 dicembre 2024 alle 04:27:19 UTC+1 Wolfgang Bangerth 
ha scritto:

> On 12/16/24 16:12, Sclah wrote:
> > 
> > I am using step-46 as starting point to implement a time dependent fsi 
> problem 
> > using block vectors and block matrices.
> > I successfully written the time iteration and changed everything to 
> block 
> > structures but I am struggling to find a method to impose the fluid 
> velocity 
> > equal to the displacement velocity at the interface (i.e. setting a fsi 
> > kinetic coupling).
> > My approach is to use AffineConstraints::add_entry and 
> > AffineConstraints::set_inhomogeneity to write something like v=(u - 
> u_old) / 
> > time_step, where u_old is the previous time step solution.
> > The problem is to find the right column index j required by these class 
> > methods. Is there a way to find it?
>
> Where do you need them, specifically? Can you use
> cell->face(f)->get_dof_indices(...)
> for example to query the DoF indices on the face of a cell, where you're 
> only 
> doing this on faces that separate the two kinds of materials?
>
> Best
> W.
>
>

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion visit 
https://groups.google.com/d/msgid/dealii/05123a90-f0ff-4e29-98cb-4467c75d66ecn%40googlegroups.com.

Reply via email to