From: Wolfgang Bangerth Date: Mon, 18 Apr 2011 23:45:22 +0000 (+0000) Subject: Minor step forward. X-Git-Tag: v8.0.0~4151 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=cf2ad1c17c2215cc5472653deeb856bcf240fda9;p=dealii.git Minor step forward. git-svn-id: https://svn.dealii.org/trunk@23606 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/examples/step-46/step-46.cc b/deal.II/examples/step-46/step-46.cc index 12194154a2..363b7cbcfc 100644 --- a/deal.II/examples/step-46/step-46.cc +++ b/deal.II/examples/step-46/step-46.cc @@ -370,9 +370,11 @@ void FluidStructureProblem::assemble_system () const unsigned int elasticity_dofs_per_cell = elasticity_fe.dofs_per_cell; FullMatrix local_matrix; + FullMatrix local_interface_matrix (elasticity_dofs_per_cell, + stokes_dofs_per_cell); Vector local_rhs; - std::vector local_dof_indices; + std::vector local_dof_indices, neighbor_dof_indices; const RightHandSide right_hand_side; @@ -461,6 +463,13 @@ void FluidStructureProblem::assemble_system () local_dof_indices.resize (dofs_per_cell); cell->get_dof_indices (local_dof_indices); + // local_rhs==0, but need to do + // this here because of + // boundary values + constraints.distribute_local_to_global (local_matrix, local_rhs, + local_dof_indices, + system_matrix, system_rhs); + // see about face terms if (cell_is_in_fluid_domain (cell)) // we are on a fluid cell @@ -484,8 +493,7 @@ void FluidStructureProblem::assemble_system () elasticity_fe_face_values.reinit (cell->neighbor(f), cell->neighbor_of_neighbor(f)); - local_matrix.reinit (stokes_dofs_per_cell, - elasticity_dofs_per_cell); + local_interface_matrix = 0; for (unsigned int q=0; q normal_vector = stokes_fe_face_values.normal_vector(q); @@ -495,13 +503,18 @@ void FluidStructureProblem::assemble_system () for (unsigned int k=0; kneighbor(f)->get_dof_indices (neighbor_dof_indices); + //constraints.distribute_local_to_global(); } } else if ((cell->neighbor(f)->level() == cell->level()) @@ -517,13 +530,6 @@ void FluidStructureProblem::assemble_system () continue; } } - // local_rhs==0, but need to do - // this here because of - // boundary values - constraints.distribute_local_to_global (local_matrix, local_rhs, - local_dof_indices, - system_matrix, system_rhs); - } }