From: Wolfgang Bangerth Date: Mon, 23 Apr 2001 15:24:10 +0000 (+0000) Subject: Implement VectorTools::interpolate_boundary_values for FunctionMap argument. X-Git-Tag: v8.0.0~19322 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f410a0b37c19c377f27cfd2dc8db41295aaeee6e;p=dealii.git Implement VectorTools::interpolate_boundary_values for FunctionMap argument. git-svn-id: https://svn.dealii.org/trunk@4464 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/numerics/vectors.h b/deal.II/deal.II/include/numerics/vectors.h index 24068ac111..0c556f34ad 100644 --- a/deal.II/deal.II/include/numerics/vectors.h +++ b/deal.II/deal.II/include/numerics/vectors.h @@ -439,24 +439,29 @@ class VectorTools const Function &rhs, Vector &rhs_vector); -//TODO:[WB] Update interpolate_boundary_values for use of FunctionMap. -// keep both, the functions using and the functions not using a FunctionMap. /** - * Prepare Dirichlet boundary conditions. - * Make up the list of nodes subject - * to Dirichlet boundary conditions - * and the values to be - * assigned to them, by interpolation around - * the boundary. If the - * @p{boundary_values} contained values - * before, the new ones are added, or - * the old one overwritten if a node - * of the boundary part to be projected - * on already was in the variable. + * Prepare Dirichlet boundary + * conditions. Make up the list + * of degrees of freedom subject + * to Dirichlet boundary + * conditions and the values to + * be assigned to them, by + * interpolation around the + * boundary. If the + * @p{boundary_values} contained + * values before, the new ones + * are added, or the old one + * overwritten if a node of the + * boundary part to be projected + * on already was in the + * variable. * - * The parameter @p{boundary_component} corresponds - * to the number @p{boundary_indicator} of the face. - * 255 is an illegal value, since it is reserved + * The parameter + * @p{boundary_component} + * corresponds to the number + * @p{boundary_indicator} of the + * face. 255 is an illegal + * value, since it is reserved * for interior faces. * * The flags in the last @@ -467,13 +472,13 @@ class VectorTools * specified by the default value * (i.e. an empty array), all * components are - * interpolated. If is different - * from the default value, it is - * assumed that the number of - * entries equals the number of - * components in the boundary - * functions and the finite - * element. + * interpolated. If it is + * different from the default + * value, it is assumed that the + * number of entries equals the + * number of components in the + * boundary functions and the + * finite element. * * It is assumed that the number * of components of the function @@ -485,15 +490,35 @@ class VectorTools * information. */ template + static void interpolate_boundary_values (const Mapping &mapping, + const DoFHandler &dof, + const typename FunctionMap::type &function_map, + std::map &boundary_values, + const std::vector &component_mask = std::vector()); + + /** + * Same function as above, but + * taking only one pair of + * boundary indicator and + * corresponding boundary + * function. Calls the other + * function with remapped + * arguments. + * + * This function is there mainly + * for backward compatibility. + */ + template static void interpolate_boundary_values (const Mapping &mapping, const DoFHandler &dof, const unsigned char boundary_component, const Function &boundary_function, std::map &boundary_values, const std::vector &component_mask = std::vector()); - + /** - * Calls the @p{interpolate_boundary_values} + * Calls the other + * @p{interpolate_boundary_values} * function, see above, with * @p{mapping=MappingQ1()}. */ @@ -504,6 +529,19 @@ class VectorTools std::map &boundary_values, const std::vector &component_mask = std::vector()); + /** + * Calls the other + * @p{interpolate_boundary_values} + * function, see above, with + * @p{mapping=MappingQ1()}. + */ + template + static void interpolate_boundary_values (const DoFHandler &dof, + const typename FunctionMap::type &function_map, + std::map &boundary_values, + const std::vector &component_mask = std::vector()); + + //TODO:[WB] Update project_boundary_values for more components /** * Project @p{function} to the boundary diff --git a/deal.II/deal.II/source/numerics/vectors.cc b/deal.II/deal.II/source/numerics/vectors.cc index 27aa6fa8e8..b21614cd81 100644 --- a/deal.II/deal.II/source/numerics/vectors.cc +++ b/deal.II/deal.II/source/numerics/vectors.cc @@ -572,7 +572,7 @@ VectorTools::interpolate_boundary_values (const Mapping<1> &, const std::vector component_mask ((component_mask_.size() == 0) ? std::vector (fe.n_components(), true) : component_mask_); - Assert (count(component_mask.begin(), component_mask.end(), true) > 0, + Assert (std::count(component_mask.begin(), component_mask.end(), true) > 0, ExcComponentMismatch()); // check whether boundary values at @@ -623,6 +623,22 @@ VectorTools::interpolate_boundary_values (const Mapping<1> &, }; + +template <> +void +VectorTools::interpolate_boundary_values (const Mapping<1> &mapping, + const DoFHandler<1> &dof, + const FunctionMap<1>::type &function_map, + std::map &boundary_values, + const std::vector &component_mask) +{ + for (FunctionMap<1>::type::const_iterator i=function_map.begin(); + i!=function_map.end(); ++i) + interpolate_boundary_values (mapping, dof, i->first, *i->second, + boundary_values, component_mask); +}; + + #endif @@ -630,19 +646,24 @@ template void VectorTools::interpolate_boundary_values (const Mapping &mapping, const DoFHandler &dof, - const unsigned char boundary_component, - const Function &boundary_function, + const typename FunctionMap::type &function_map, std::map &boundary_values, const std::vector &component_mask_) { - Assert (boundary_component != 255, + // if for whatever we were passed + // an empty map, return immediately + if (function_map.size() == 0) + return; + + Assert (function_map.find(255) == function_map.end(), ExcInvalidBoundaryIndicator()); const FiniteElement &fe = dof.get_fe(); const unsigned int n_components = fe.n_components(); const bool fe_is_system = (n_components != 1); - Assert (n_components == boundary_function.n_components, + Assert ((function_map.size() == 0) || + (n_components == function_map.begin()->second->n_components), ExcInvalidFE()); // set the component mask to either @@ -651,7 +672,7 @@ VectorTools::interpolate_boundary_values (const Mapping &mapping const std::vector component_mask ((component_mask_.size() == 0) ? std::vector (fe.n_components(), true) : component_mask_); - Assert (count(component_mask.begin(), component_mask.end(), true) > 0, + Assert (std::count(component_mask.begin(), component_mask.end(), true) > 0, ExcComponentMismatch()); // field to store the indices @@ -686,14 +707,15 @@ VectorTools::interpolate_boundary_values (const Mapping &mapping Quadrature aux_quad (unit_support_points, dummy_weights); FEFaceValues fe_values (mapping, fe, aux_quad, update_q_points); - DoFHandler::active_cell_iterator cell = dof.begin_active(), - endc = dof.end(); + typename DoFHandler::active_cell_iterator cell = dof.begin_active(), + endc = dof.end(); for (; cell!=endc; ++cell) for (unsigned int face_no = 0; face_no < GeometryInfo::faces_per_cell; ++face_no) { - DoFHandler::face_iterator face = cell->face(face_no); - if (boundary_component == face->boundary_indicator()) + typename DoFHandler::face_iterator face = cell->face(face_no); + const unsigned char boundary_component = face->boundary_indicator(); + if (function_map.find(boundary_component) != function_map.end()) // face is of the right component { // get indices, physical location and @@ -705,10 +727,10 @@ VectorTools::interpolate_boundary_values (const Mapping &mapping if (fe_is_system) { - boundary_function.vector_value_list (dof_locations, dof_values_system); + function_map.find(boundary_component)->second->vector_value_list (dof_locations, + dof_values_system); // enter into list - for (unsigned int i=0; i &mapping { // get only the one component that // this function has - boundary_function.value_list (dof_locations, - dof_values_scalar, - 0); + function_map.find(boundary_component)->second->value_list (dof_locations, + dof_values_scalar, + 0); // enter into list @@ -733,6 +755,24 @@ VectorTools::interpolate_boundary_values (const Mapping &mapping } } + + +template +void +VectorTools::interpolate_boundary_values (const Mapping &mapping, + const DoFHandler &dof, + const unsigned char boundary_component, + const Function &boundary_function, + std::map &boundary_values, + const std::vector &component_mask) +{ + typename FunctionMap::type function_map; + function_map[boundary_component] = &boundary_function; + interpolate_boundary_values (mapping, dof, function_map, boundary_values, + component_mask); +}; + + template void @@ -750,6 +790,21 @@ VectorTools::interpolate_boundary_values (const DoFHandler &dof, +template +void +VectorTools::interpolate_boundary_values (const DoFHandler &dof, + const typename FunctionMap::type &function_map, + std::map &boundary_values, + const std::vector &component_mask) +{ + Assert (DEAL_II_COMPAT_MAPPING, ExcCompatibility("mapping")); + static const MappingQ1 mapping; + interpolate_boundary_values(mapping, dof, function_map, + boundary_values, component_mask); +} + + + template void VectorTools::project_boundary_values (const Mapping &mapping, diff --git a/deal.II/doc/news/2001/c-3-1.html b/deal.II/doc/news/2001/c-3-1.html index 6fefb46dc0..169aeb6aa5 100644 --- a/deal.II/doc/news/2001/c-3-1.html +++ b/deal.II/doc/news/2001/c-3-1.html @@ -313,6 +313,16 @@ documentation, etc.

deal.II

    +
  1. + New: Implement the VectorTools::interpolate_boundary_values + a second time, this time taking a FunctionMap object as do all the other + functions of this type. +
    + (WB 2001/04/23) +

    +
  2. Fixed: The vertex numbers written by the GridOut::