Dear All,

I have a class which represents a velocity field (I am using a variant of 
the pressure correction navier stokes model with two phases) and a class 
which represents a transport model along a provided velocity field.   I am 
adding in some terms into the velocity model (currently surface tension but 
soon to be added is the canham-helfrich or willmore energy) accounting for 
the physics and force balance on the boundary between the two phases.  The 
approach that we are using is to create a block structure inside the 
velocity class to incorporate these other physics into the model.  We have 
a splitting of velocity separate from the pressure terms and have separate 
dofhandlers for these two groups but now we are lumping all the components 
that need to be simultaneously solved for with the velocity into the 
velocity dofhandler and fesystem.  For instance we have the block 
components in the velocity fesystem:

dim components of velocity,  normal velocity, curvature and willmore energy

so that our velocity fesystem looks like

FESystem<dim> (

                             FESystem<dim>(FE_Q<dim> 
(parameters.degree_of_fe_velocity),dim), 1,  // velocity

                             FE_Q<dim>(parameters.degree_of_fe_velocity), 
1, // normal velocity

                             FE_Q<dim>(parameters.degree_of_fe_velocity), 
1, // curvature

                             FE_Q<dim>(parameters.degree_of_fe_velocity), 1 
// willmore force

                             )

I have figured out how to construct the block sparsity patterns and 
matrices etc and solve for this, but at the end of the day this class 
represents a velocity not a velocity plus other stuff.  This class 
interfaces with the other class representing transport and I need to be 
able to provide a velocity, dofhandler and fesystem in such a way that I 
don't break encapsulation ie the other class doesn't need to know that 
internally there was a block system and other components... 

The other class expects to be able to query a 
TrilinosWrappers::MPI::Vector,  FESystem<dim> and DoFHandler<dim> 
 representing the velocity finite element solution.

Supposing internally to the velocity class I have

solution, a TrilinosWrapper::MPI::BlockVector
fe_system, the FESystem<dim> described above and 
dof_handler, the associated DoFHandler<dim> which is sorted 
DoFRenumbering::hierarchical and DoFRenumbering::component_wise.

Supposing velocity_block represents the block number associated to the 
velocity components, I noticed that solution.block{velocity_block} would be 
a TrilinosWrapper::MPI::Vector of the velocity as desired and it is 
possible to extract the sub FiniteElement 
fe_system.base_element{velocity_block} corresponding only to the velocity, 
but the dof_handler has no such existing capability and if I pass it as is 
with the extracted FiniteElement into an FEValues, an assertion is 
triggered when we call fe_values.reinit(cell) saying that the given finite 
element does not match the dofhandler's finite element...  I don't want to 
pass the entire system since that requires the other class to know that it 
only needs the first block components etc.

My question is this:  What would be a best approach to dealing with this? 
 Is there currently a way to extract a "diagonal block" from the system and 
reduce the solution, fesystem and dofhandler to only represent that part? 
 Do I need to create a new FESystem, DoFHandler and vector to provide the 
functionality I need ?

I hope this is clear enough and enough(and not too much) background to 
understand what I am asking.  Thanks for any help!

Spencer Patty





-- 
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.
For more options, visit https://groups.google.com/d/optout.

Reply via email to