From 97cb6e17e19f308fe9bf70b853047271e8c1fd83 Mon Sep 17 00:00:00 2001 From: Ralf Hartmann Date: Tue, 18 Apr 2000 13:20:43 +0000 Subject: [PATCH] remove Assembler from create_right_hand_side function. New code for that function for scalar and system FEs. git-svn-id: https://svn.dealii.org/trunk@2744 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/deal.II/source/numerics/vectors.cc | 90 ++++++++++++++++++---- 1 file changed, 73 insertions(+), 17 deletions(-) diff --git a/deal.II/deal.II/source/numerics/vectors.cc b/deal.II/deal.II/source/numerics/vectors.cc index e988a14fbf..c2587cb979 100644 --- a/deal.II/deal.II/source/numerics/vectors.cc +++ b/deal.II/deal.II/source/numerics/vectors.cc @@ -456,33 +456,89 @@ void VectorTools::project (const DoFHandler &dof, template -void VectorTools::create_right_hand_side (const DoFHandler &dof, +void VectorTools::create_right_hand_side (const DoFHandler &dof_handler, const Quadrature &quadrature, - const Function &rhs, + const Function &rhs_function, Vector &rhs_vector) { - Assert (dof.get_fe().n_components() == rhs.n_components, + const FiniteElement &fe = dof_handler.get_fe(); + Assert (fe.n_components() == rhs_function.n_components, ExcComponentMismatch()); + Assert (rhs_vector.size() == dof_handler.n_dofs(), + ExcDimensionMismatch(rhs_vector.size(), dof_handler.n_dofs())); + rhs_vector.clear(); UpdateFlags update_flags = UpdateFlags(update_values | update_q_points | update_JxW_values); - SparseMatrix dummy; - const Assembler::AssemblerData data (dof, - false, true, - dummy, rhs_vector, - quadrature, update_flags); - TriaActiveIterator > - assembler (const_cast*>(&dof.get_tria()), - dof.get_tria().begin_active()->level(), - dof.get_tria().begin_active()->index(), - &data); - MassMatrix equation(&rhs,0); - do + FEValues fe_values (fe, quadrature, update_flags); + + const unsigned int dofs_per_cell = fe_values.dofs_per_cell, + n_q_points = fe_values.n_quadrature_points, + n_components = fe.n_components(); + + vector dofs (dofs_per_cell); + Vector cell_vector (dofs_per_cell); + + DoFHandler::active_cell_iterator cell = dof_handler.begin_active(), + endc = dof_handler.end(); + + if (n_components==1) { - assembler->assemble (equation); + vector rhs_values(n_q_points); + + for (; cell!=endc; ++cell) + { + fe_values.reinit(cell); + + const FullMatrix &values = fe_values.get_shape_values (); + const vector &weights = fe_values.get_JxW_values (); + rhs_function.value_list (fe_values.get_quadrature_points(), rhs_values); + + cell_vector.clear(); + for (unsigned int point=0; pointget_dof_indices (dofs); + + for (unsigned int i=0; i > rhs_values(n_q_points, Vector(n_components)); + + for (; cell!=endc; ++cell) + { + fe_values.reinit(cell); + + const FullMatrix &values = fe_values.get_shape_values (); + const vector &weights = fe_values.get_JxW_values (); + rhs_function.vector_value_list (fe_values.get_quadrature_points(), rhs_values); + + cell_vector.clear(); + for (unsigned int point=0; pointget_dof_indices (dofs); + + for (unsigned int i=0; i