From: Daniel Garcia-Sanchez Date: Sat, 22 Jun 2019 08:44:47 +0000 (+0200) Subject: Add std::complex project_boundary_values_curl_conforming_l2 X-Git-Tag: v9.2.0-rc1~1425^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f11f0a7015255e76eb5e258f0565ba7de3db0fcb;p=dealii.git Add std::complex project_boundary_values_curl_conforming_l2 --- diff --git a/include/deal.II/numerics/vector_tools.h b/include/deal.II/numerics/vector_tools.h index 1094d36f72..bff31509d9 100644 --- a/include/deal.II/numerics/vector_tools.h +++ b/include/deal.II/numerics/vector_tools.h @@ -1682,14 +1682,14 @@ namespace VectorTools * @see * @ref GlossBoundaryIndicator "Glossary entry on boundary indicators" */ - template + template void project_boundary_values_curl_conforming_l2( const DoFHandler & dof_handler, const unsigned int first_vector_component, - const Function &boundary_function, + const Function &boundary_function, const types::boundary_id boundary_component, - AffineConstraints & constraints, + AffineConstraints & constraints, const Mapping & mapping = StaticMappingQ1::mapping); @@ -1699,14 +1699,14 @@ namespace VectorTools * * @ingroup constraints */ - template + template void project_boundary_values_curl_conforming_l2( const hp::DoFHandler & dof_handler, const unsigned int first_vector_component, - const Function & boundary_function, + const Function & boundary_function, const types::boundary_id boundary_component, - AffineConstraints & constraints, + AffineConstraints & constraints, const hp::MappingCollection &mapping_collection = hp::StaticMappingQ1::mapping_collection); diff --git a/include/deal.II/numerics/vector_tools.templates.h b/include/deal.II/numerics/vector_tools.templates.h index 2798737c53..0e49213734 100644 --- a/include/deal.II/numerics/vector_tools.templates.h +++ b/include/deal.II/numerics/vector_tools.templates.h @@ -4374,16 +4374,16 @@ namespace VectorTools // projection of the boundary // function on the interior of // faces. - template + template void compute_face_projection_curl_conforming( - const cell_iterator &cell, - const unsigned int face, - hp::FEValues & hp_fe_values, - const Function &boundary_function, - const unsigned int first_vector_component, - std::vector &dof_values, - std::vector & dofs_processed) + const cell_iterator & cell, + const unsigned int face, + hp::FEValues & hp_fe_values, + const Function &boundary_function, + const unsigned int first_vector_component, + std::vector & dof_values, + std::vector & dofs_processed) { const unsigned int spacedim = dim; hp_fe_values.reinit(cell, @@ -5315,15 +5315,15 @@ namespace VectorTools namespace internals { - template - void - compute_edge_projection_l2(const cell_iterator &cell, - const unsigned int face, - const unsigned int line, - hp::FEValues<3> & hp_fe_values, - const Function<3> & boundary_function, + template + typename std::enable_if::type + compute_edge_projection_l2(const cell_iterator & cell, + const unsigned int face, + const unsigned int line, + hp::FEValues & hp_fe_values, + const Function &boundary_function, const unsigned int first_vector_component, - std::vector &dof_values, + std::vector &dof_values, std::vector & dofs_processed) { // This function computes the L2-projection of the given @@ -5333,8 +5333,7 @@ namespace VectorTools // In the context of this function, by associated DoFs we mean: // the DoFs corresponding to the group of components making up the vector // with first component first_vector_component (length dim). - const unsigned int dim = 3; - const FiniteElement &fe = cell->get_fe(); + const FiniteElement &fe = cell->get_fe(); // reinit for this cell, face and line. hp_fe_values.reinit( @@ -5351,8 +5350,8 @@ namespace VectorTools const std::vector> &quadrature_points = fe_values.get_quadrature_points(); - std::vector> values(fe_values.n_quadrature_points, - Vector(fe.n_components())); + std::vector> values(fe_values.n_quadrature_points, + Vector(fe.n_components())); // Get boundary function values // at quadrature points. @@ -5456,10 +5455,10 @@ namespace VectorTools // Matrix and RHS vectors to store linear system: // We have (degree+1) basis functions for an edge - FullMatrix edge_matrix(degree + 1, degree + 1); - FullMatrix edge_matrix_inv(degree + 1, degree + 1); - Vector edge_rhs(degree + 1); - Vector edge_solution(degree + 1); + FullMatrix edge_matrix(degree + 1, degree + 1); + FullMatrix edge_matrix_inv(degree + 1, degree + 1); + Vector edge_rhs(degree + 1); + Vector edge_solution(degree + 1); const FEValuesExtractors::Vector vec(first_vector_component); @@ -5560,15 +5559,15 @@ namespace VectorTools } - template - void + template + typename std::enable_if::type compute_edge_projection_l2(const cell_iterator &, const unsigned int, const unsigned int, hp::FEValues &, - const Function &, + const Function &, const unsigned int, - std::vector &, + std::vector &, std::vector &) { // dummy implementation of above function @@ -5576,16 +5575,16 @@ namespace VectorTools Assert(false, ExcInternalError()); } - template + template void compute_face_projection_curl_conforming_l2( - const cell_iterator & cell, - const unsigned int face, - hp::FEFaceValues &hp_fe_face_values, - const Function & boundary_function, - const unsigned int first_vector_component, - std::vector & dof_values, - std::vector & dofs_processed) + const cell_iterator & cell, + const unsigned int face, + hp::FEFaceValues & hp_fe_face_values, + const Function &boundary_function, + const unsigned int first_vector_component, + std::vector & dof_values, + std::vector & dofs_processed) { // This function computes the L2-projection of the boundary // function on the interior of faces only. In 3D, this should only be @@ -5608,8 +5607,8 @@ namespace VectorTools fe_face_values.get_quadrature_points(); const unsigned int degree = fe.degree - 1; - std::vector> values(fe_face_values.n_quadrature_points, - Vector(fe.n_components())); + std::vector> values(fe_face_values.n_quadrature_points, + Vector(fe.n_components())); // Get boundary function values at quadrature points. boundary_function.vector_value_list(quadrature_points, values); @@ -5690,10 +5689,10 @@ namespace VectorTools // Matrix and RHS vectors to store: // We have (degree+1) edge basis functions - FullMatrix edge_matrix(degree + 1, degree + 1); - FullMatrix edge_matrix_inv(degree + 1, degree + 1); - Vector edge_rhs(degree + 1); - Vector edge_solution(degree + 1); + FullMatrix edge_matrix(degree + 1, degree + 1); + FullMatrix edge_matrix_inv(degree + 1, degree + 1); + Vector edge_rhs(degree + 1); + Vector edge_solution(degree + 1); const FEValuesExtractors::Vector vec(first_vector_component); @@ -5893,19 +5892,19 @@ namespace VectorTools // There are 2*degree*(degree+1) DoFs associated with a face in // 3D. Note this doesn't include the DoFs associated with edges on // that face. - FullMatrix face_matrix(2 * degree * (degree + 1)); - FullMatrix face_matrix_inv(2 * degree * (degree + 1)); - Vector face_rhs(2 * degree * (degree + 1)); - Vector face_solution(2 * degree * (degree + 1)); + FullMatrix face_matrix(2 * degree * (degree + 1)); + FullMatrix face_matrix_inv(2 * degree * (degree + 1)); + Vector face_rhs(2 * degree * (degree + 1)); + Vector face_solution(2 * degree * (degree + 1)); // Project the boundary function onto the shape functions for this // face and set up a linear system of equations to get the values // for the DoFs associated with this face. We also must include // the residuals from the shape functions associated with edges. - Tensor<1, dim> tmp; - Tensor<1, dim> cross_product_i; - Tensor<1, dim> cross_product_j; - Tensor<1, dim> cross_product_rhs; + Tensor<1, dim, number> tmp; + Tensor<1, dim> cross_product_i; + Tensor<1, dim> cross_product_j; + Tensor<1, dim, number> cross_product_rhs; // Loop to construct face linear system. for (unsigned int q_point = 0; @@ -6017,14 +6016,14 @@ namespace VectorTools } } - template + template void compute_project_boundary_values_curl_conforming_l2( const DoFHandlerType & dof_handler, const unsigned int first_vector_component, - const Function & boundary_function, + const Function & boundary_function, const types::boundary_id boundary_component, - AffineConstraints & constraints, + AffineConstraints & constraints, const hp::MappingCollection &mapping_collection) { // L2-projection based interpolation formed in one (in 2D) or two (in 3D) @@ -6079,7 +6078,7 @@ namespace VectorTools // Storage for dof values found and whether they have been processed: std::vector dofs_processed; - std::vector dof_values; + std::vector dof_values; std::vector face_dof_indices; typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active(); @@ -6332,15 +6331,15 @@ namespace VectorTools } // namespace internals - template + template void project_boundary_values_curl_conforming_l2( - const DoFHandler & dof_handler, - const unsigned int first_vector_component, - const Function & boundary_function, - const types::boundary_id boundary_component, - AffineConstraints &constraints, - const Mapping & mapping) + const DoFHandler & dof_handler, + const unsigned int first_vector_component, + const Function &boundary_function, + const types::boundary_id boundary_component, + AffineConstraints & constraints, + const Mapping & mapping) { // non-hp version - calls the internal // compute_project_boundary_values_curl_conforming_l2() function @@ -6356,14 +6355,14 @@ namespace VectorTools mapping_collection); } - template + template void project_boundary_values_curl_conforming_l2( const hp::DoFHandler & dof_handler, const unsigned int first_vector_component, - const Function & boundary_function, + const Function & boundary_function, const types::boundary_id boundary_component, - AffineConstraints & constraints, + AffineConstraints & constraints, const hp::MappingCollection &mapping_collection) { // hp version - calls the internal diff --git a/source/numerics/vector_tools_boundary.inst.in b/source/numerics/vector_tools_boundary.inst.in index c4f88aa271..0486d98d57 100644 --- a/source/numerics/vector_tools_boundary.inst.in +++ b/source/numerics/vector_tools_boundary.inst.in @@ -220,22 +220,42 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) AffineConstraints &, const hp::MappingCollection &); template void - project_boundary_values_curl_conforming_l2( + project_boundary_values_curl_conforming_l2( const DoFHandler &, const unsigned int, - const Function &, + const Function &, const types::boundary_id, AffineConstraints &, const Mapping &); +# ifdef DEAL_II_WITH_COMPLEX_VALUES + template void + project_boundary_values_curl_conforming_l2( + const DoFHandler &, + const unsigned int, + const Function> &, + const types::boundary_id, + AffineConstraints> &, + const Mapping &); +# endif template void - project_boundary_values_curl_conforming_l2( + project_boundary_values_curl_conforming_l2( const hp::DoFHandler &, const unsigned int, - const Function &, + const Function &, const types::boundary_id, AffineConstraints &, const hp::MappingCollection &); template void +# ifdef DEAL_II_WITH_COMPLEX_VALUES + project_boundary_values_curl_conforming_l2( + const hp::DoFHandler &, + const unsigned int, + const Function> &, + const types::boundary_id, + AffineConstraints> &, + const hp::MappingCollection &); +# endif + template void project_boundary_values_div_conforming( const DoFHandler &, const unsigned int,