From: Wolfgang Bangerth Date: Thu, 22 Apr 2021 19:35:01 +0000 (-0600) Subject: Simplify some type computations in FEValuesViews. X-Git-Tag: v9.3.0-rc1~179^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=a52daa9f71fb875df4ebc082ac0f2f42a3fc332a;p=dealii.git Simplify some type computations in FEValuesViews. Specifically, instead of using a nested class, just use a template type alias. While there, also fix a few places where we should have used the existing type alias but instead had repeated the same type computation used in the definition of the type alias. --- diff --git a/include/deal.II/fe/fe_values.h b/include/deal.II/fe/fe_values.h index 8b5a704d9a..910717389a 100644 --- a/include/deal.II/fe/fe_values.h +++ b/include/deal.II/fe/fe_values.h @@ -175,52 +175,53 @@ namespace FEValuesViews using third_derivative_type = dealii::Tensor<3, spacedim>; /** - * A struct that provides the output type for the product of the value - * and derivatives of basis functions of the Scalar view and any @p Number type. + * An alias for the data type of the product of a @p Number and the + * values of the view this class provides. This is the data type of + * scalar components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. */ template - struct OutputType - { - /** - * An alias for the data type of the product of a @p Number and the - * values of the view the Scalar class. - */ - using value_type = - typename ProductType::value_type>::type; + using solution_value_type = typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * gradients of the view the Scalar class. - */ - using gradient_type = typename ProductType< - Number, - typename Scalar::gradient_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * gradients of the view this class provides. This is the data type of + * scalar components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_gradient_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * laplacians of the view the Scalar class. - */ - using laplacian_type = - typename ProductType::value_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * laplacians of the view this class provides. This is the data type of + * scalar components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_laplacian_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * hessians of the view the Scalar class. - */ - using hessian_type = typename ProductType< - Number, - typename Scalar::hessian_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * hessians of the view this class provides. This is the data type of + * scalar components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_hessian_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * third derivatives of the view the Scalar class. - */ - using third_derivative_type = typename ProductType< - Number, - typename Scalar::third_derivative_type>::type; - }; + /** + * An alias for the data type of the product of a @p Number and the + * third derivatives of the view this class provides. This is the data type + * of scalar components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_third_derivative_type = + typename ProductType::type; /** * A structure where for each shape function we pre-compute a bunch of @@ -348,8 +349,7 @@ namespace FEValuesViews void get_function_values( const InputVector &fe_function, - std::vector::type> + std::vector> &values) const; /** @@ -390,8 +390,7 @@ namespace FEValuesViews void get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> + std::vector> &values) const; /** @@ -415,8 +414,7 @@ namespace FEValuesViews void get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const; /** @@ -429,8 +427,7 @@ namespace FEValuesViews void get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const; /** @@ -454,8 +451,7 @@ namespace FEValuesViews void get_function_hessians( const InputVector &fe_function, - std::vector::type> + std::vector> &hessians) const; /** @@ -468,8 +464,7 @@ namespace FEValuesViews void get_function_hessians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::hessian_type> + std::vector> &hessians) const; @@ -495,8 +490,7 @@ namespace FEValuesViews void get_function_laplacians( const InputVector &fe_function, - std::vector::type> + std::vector> &laplacians) const; /** @@ -509,8 +503,7 @@ namespace FEValuesViews void get_function_laplacians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::laplacian_type> + std::vector> &laplacians) const; @@ -536,8 +529,8 @@ namespace FEValuesViews void get_function_third_derivatives( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_third_derivative_type> &third_derivatives) const; /** @@ -549,9 +542,10 @@ namespace FEValuesViews template void get_function_third_derivatives_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - third_derivative_type> &third_derivatives) const; + const InputVector &dof_values, + std::vector< + solution_third_derivative_type> + &third_derivatives) const; private: @@ -666,76 +660,82 @@ namespace FEValuesViews using third_derivative_type = dealii::Tensor<4, spacedim>; /** - * A struct that provides the output type for the product of the value - * and derivatives of basis functions of the Vector view and any @p Number type. + * An alias for the data type of the product of a @p Number and the + * values of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. */ template - struct OutputType - { - /** - * An alias for the data type of the product of a @p Number and the - * values of the view the Vector class. - */ - using value_type = - typename ProductType::value_type>::type; + using solution_value_type = typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * gradients of the view the Vector class. - */ - using gradient_type = typename ProductType< - Number, - typename Vector::gradient_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * gradients of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_gradient_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * symmetric gradients of the view the Vector class. - */ - using symmetric_gradient_type = typename ProductType< - Number, - typename Vector::symmetric_gradient_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * symmetric gradients of the view this class provides. This is the data + * type of vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_symmetric_gradient_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * divergences of the view the Vector class. - */ - using divergence_type = typename ProductType< - Number, - typename Vector::divergence_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * divergences of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_divergence_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * laplacians of the view the Vector class. - */ - using laplacian_type = - typename ProductType::value_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * laplacians of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_laplacian_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * curls of the view the Vector class. - */ - using curl_type = - typename ProductType::curl_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * curls of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_curl_type = typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * hessians of the view the Vector class. - */ - using hessian_type = typename ProductType< - Number, - typename Vector::hessian_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * hessians of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_hessian_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * third derivatives of the view the Vector class. - */ - using third_derivative_type = typename ProductType< - Number, - typename Vector::third_derivative_type>::type; - }; + /** + * An alias for the data type of the product of a @p Number and the + * third derivatives of the view this class provides. This is the data type + * of vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_third_derivative_type = + typename ProductType::type; /** * A structure where for each shape function we pre-compute a bunch of @@ -940,8 +940,7 @@ namespace FEValuesViews void get_function_values( const InputVector &fe_function, - std::vector::type> + std::vector> &values) const; /** @@ -982,8 +981,7 @@ namespace FEValuesViews void get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> + std::vector> &values) const; /** @@ -1007,8 +1005,7 @@ namespace FEValuesViews void get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const; /** @@ -1021,8 +1018,7 @@ namespace FEValuesViews void get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const; /** @@ -1052,8 +1048,8 @@ namespace FEValuesViews void get_function_symmetric_gradients( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_symmetric_gradient_type> &symmetric_gradients) const; /** @@ -1065,9 +1061,10 @@ namespace FEValuesViews template void get_function_symmetric_gradients_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - symmetric_gradient_type> &symmetric_gradients) const; + const InputVector &dof_values, + std::vector< + solution_symmetric_gradient_type> + &symmetric_gradients) const; /** * Return the divergence of the selected vector components of the finite @@ -1091,8 +1088,7 @@ namespace FEValuesViews void get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const; /** @@ -1105,8 +1101,7 @@ namespace FEValuesViews void get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const; /** @@ -1131,9 +1126,8 @@ namespace FEValuesViews void get_function_curls( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &curls) const; + std::vector> &curls) + const; /** * This function relates to get_function_curls() in the same way @@ -1145,9 +1139,8 @@ namespace FEValuesViews void get_function_curls_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::curl_type> - &curls) const; + std::vector> &curls) + const; /** * Return the Hessians of the selected vector components of the finite @@ -1170,8 +1163,7 @@ namespace FEValuesViews void get_function_hessians( const InputVector &fe_function, - std::vector::type> + std::vector> &hessians) const; /** @@ -1184,8 +1176,7 @@ namespace FEValuesViews void get_function_hessians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::hessian_type> + std::vector> &hessians) const; /** @@ -1210,8 +1201,7 @@ namespace FEValuesViews void get_function_laplacians( const InputVector &fe_function, - std::vector::type> + std::vector> &laplacians) const; /** @@ -1224,8 +1214,7 @@ namespace FEValuesViews void get_function_laplacians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::laplacian_type> + std::vector> &laplacians) const; /** @@ -1250,8 +1239,8 @@ namespace FEValuesViews void get_function_third_derivatives( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_third_derivative_type> &third_derivatives) const; /** @@ -1263,9 +1252,10 @@ namespace FEValuesViews template void get_function_third_derivatives_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - third_derivative_type> &third_derivatives) const; + const InputVector &dof_values, + std::vector< + solution_third_derivative_type> + &third_derivatives) const; private: /** @@ -1335,28 +1325,23 @@ namespace FEValuesViews using divergence_type = dealii::Tensor<1, spacedim>; /** - * A struct that provides the output type for the product of the value - * and derivatives of basis functions of the SymmetricTensor view and any @p Number type. + * An alias for the data type of the product of a @p Number and the + * values of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. */ template - struct OutputType - { - /** - * An alias for the data type of the product of a @p Number and the - * values of the view the SymmetricTensor class. - */ - using value_type = typename ProductType< - Number, - typename SymmetricTensor<2, dim, spacedim>::value_type>::type; + using solution_value_type = typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * divergences of the view the SymmetricTensor class. - */ - using divergence_type = typename ProductType< - Number, - typename SymmetricTensor<2, dim, spacedim>::divergence_type>::type; - }; + /** + * An alias for the data type of the product of a @p Number and the + * divergences of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_divergence_type = + typename ProductType::type; /** * A structure where for each shape function we pre-compute a bunch of @@ -1484,8 +1469,7 @@ namespace FEValuesViews void get_function_values( const InputVector &fe_function, - std::vector::type> + std::vector> &values) const; /** @@ -1526,8 +1510,7 @@ namespace FEValuesViews void get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> + std::vector> &values) const; /** @@ -1555,8 +1538,7 @@ namespace FEValuesViews void get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const; /** @@ -1569,8 +1551,7 @@ namespace FEValuesViews void get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const; private: @@ -1635,36 +1616,33 @@ namespace FEValuesViews using gradient_type = dealii::Tensor<3, spacedim>; /** - * A struct that provides the output type for the product of the value - * and derivatives of basis functions of the Tensor view and any @p Number type. + * An alias for the data type of the product of a @p Number and the + * values of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. */ template - struct OutputType - { - /** - * An alias for the data type of the product of a @p Number and the - * values of the view the Tensor class. - */ - using value_type = typename ProductType< - Number, - typename Tensor<2, dim, spacedim>::value_type>::type; + using solution_value_type = typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * divergences of the view the Tensor class. - */ - using divergence_type = typename ProductType< - Number, - typename Tensor<2, dim, spacedim>::divergence_type>::type; + /** + * An alias for the data type of the product of a @p Number and the + * divergences of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_divergence_type = + typename ProductType::type; - /** - * An alias for the data type of the product of a @p Number and the - * gradient of the view the Tensor class. - */ - using gradient_type = typename ProductType< - Number, - typename Tensor<2, dim, spacedim>::gradient_type>::type; - }; + /** + * An alias for the data type of the product of a @p Number and the + * gradient of the view this class provides. This is the data type of + * vector components of a finite element field whose degrees of + * freedom are described by a vector with elements of type @p Number. + */ + template + using solution_gradient_type = + typename ProductType::type; /** * A structure where for each shape function we pre-compute a bunch of @@ -1809,8 +1787,7 @@ namespace FEValuesViews void get_function_values( const InputVector &fe_function, - std::vector::type> + std::vector> &values) const; /** @@ -1851,8 +1828,7 @@ namespace FEValuesViews void get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> + std::vector> &values) const; /** @@ -1880,8 +1856,7 @@ namespace FEValuesViews void get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const; /** @@ -1894,8 +1869,7 @@ namespace FEValuesViews void get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const; /** @@ -1918,8 +1892,7 @@ namespace FEValuesViews void get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const; /** @@ -1932,8 +1905,7 @@ namespace FEValuesViews void get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const; private: diff --git a/include/deal.II/meshworker/scratch_data.h b/include/deal.II/meshworker/scratch_data.h index 36ab1e1088..9f63419326 100644 --- a/include/deal.II/meshworker/scratch_data.h +++ b/include/deal.II/meshworker/scratch_data.h @@ -36,17 +36,6 @@ DEAL_II_NAMESPACE_OPEN namespace MeshWorker { - namespace internal - { - /** - * Alias used to deduce the OutputType when asking for values or gradients - * of shape functions multiplied with the given @p NumberType. - */ - template - using OutputType = typename FEValuesViews::View:: - template OutputType; - } // namespace internal - /** * A helper class to simplify the parallel assembly of linear and non-linear * problems, and the evaluation of finite element fields. @@ -606,9 +595,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - value_type> & + const std::vector:: + template solution_value_type> & get_values(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -635,9 +623,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - gradient_type> & + const std::vector:: + template solution_gradient_type> & get_gradients(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -665,9 +652,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - symmetric_gradient_type> & + const std::vector:: + template solution_symmetric_gradient_type> & get_symmetric_gradients(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -694,9 +680,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - divergence_type> & + const std::vector:: + template solution_divergence_type> & get_divergences(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -723,12 +708,11 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector::curl_type> - & - get_curls(const std::string &global_vector_name, - const Extractor & variable, - const Number dummy = Number(0)); + const std::vector:: + template solution_curl_type> & + get_curls(const std::string &global_vector_name, + const Extractor & variable, + const Number dummy = Number(0)); /** * For the solution vector identified by @p global_vector_name, compute @@ -752,9 +736,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - hessian_type> & + const std::vector:: + template solution_hessian_type> & get_hessians(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -781,9 +764,8 @@ namespace MeshWorker * the documentation of the FEValues class for more information. */ template - const std::vector< - typename internal::OutputType:: - laplacian_type> & + const std::vector:: + template solution_laplacian_type> & get_laplacians(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -810,9 +792,8 @@ namespace MeshWorker * FiniteElement" in the documentation of the FEValues for more information. */ template - const std::vector< - typename internal::OutputType:: - third_derivative_type> & + const std::vector:: + template solution_third_derivative_type> & get_third_derivatives(const std::string &global_vector_name, const Extractor & variable, const Number dummy = Number(0)); @@ -1054,13 +1035,11 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType::value_type> - & - ScratchData::get_values( - const std::string &global_vector_name, - const Extractor & variable, - const Number dummy) + const std::vector:: + template solution_value_type> & + ScratchData::get_values(const std::string &global_vector_name, + const Extractor & variable, + const Number dummy) { const std::vector &independent_local_dofs = get_local_dof_values(global_vector_name, dummy); @@ -1074,8 +1053,8 @@ namespace MeshWorker // Now build the return type using RetType = - std::vector::value_type>; + std::vector:: + template solution_value_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1092,9 +1071,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - gradient_type> & + const std::vector:: + template solution_gradient_type> & ScratchData::get_gradients( const std::string &global_vector_name, const Extractor & variable, @@ -1111,9 +1089,9 @@ namespace MeshWorker global_vector_name, variable, "_gradients_q", n_q_points, dummy); // Now build the return type - using RetType = std::vector< - typename internal::OutputType:: - gradient_type>; + using RetType = + std::vector:: + template solution_gradient_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1130,9 +1108,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - hessian_type> & + const std::vector:: + template solution_hessian_type> & ScratchData::get_hessians( const std::string &global_vector_name, const Extractor & variable, @@ -1150,8 +1127,8 @@ namespace MeshWorker // Now build the return type using RetType = - std::vector::hessian_type>; + std::vector:: + template solution_hessian_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1169,9 +1146,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - laplacian_type> & + const std::vector:: + template solution_laplacian_type> & ScratchData::get_laplacians( const std::string &global_vector_name, const Extractor & variable, @@ -1188,9 +1164,9 @@ namespace MeshWorker global_vector_name, variable, "_laplacians_q", n_q_points, dummy); // Now build the return type - using RetType = std::vector< - typename internal::OutputType:: - laplacian_type>; + using RetType = + std::vector:: + template solution_laplacian_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1208,9 +1184,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - third_derivative_type> & + const std::vector:: + template solution_third_derivative_type> & ScratchData::get_third_derivatives( const std::string &global_vector_name, const Extractor & variable, @@ -1227,9 +1202,9 @@ namespace MeshWorker global_vector_name, variable, "_third_derivatives_q", n_q_points, dummy); // Now build the return type - using RetType = std::vector< - typename internal::OutputType:: - third_derivative_type>; + using RetType = + std::vector:: + template solution_third_derivative_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1247,9 +1222,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - symmetric_gradient_type> & + const std::vector:: + template solution_symmetric_gradient_type> & ScratchData::get_symmetric_gradients( const std::string &global_vector_name, const Extractor & variable, @@ -1267,9 +1241,9 @@ namespace MeshWorker // Now build the return type - using RetType = std::vector< - typename internal::OutputType:: - symmetric_gradient_type>; + using RetType = + std::vector:: + template solution_symmetric_gradient_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1286,9 +1260,8 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType:: - divergence_type> & + const std::vector:: + template solution_divergence_type> & ScratchData::get_divergences( const std::string &global_vector_name, const Extractor & variable, @@ -1305,9 +1278,9 @@ namespace MeshWorker global_vector_name, variable, "_divergence_q", n_q_points, dummy); // Now build the return type - using RetType = std::vector< - typename internal::OutputType:: - divergence_type>; + using RetType = + std::vector:: + template solution_divergence_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( @@ -1325,12 +1298,11 @@ namespace MeshWorker template template - const std::vector< - typename internal::OutputType::curl_type> - & - ScratchData::get_curls(const std::string &global_vector_name, - const Extractor & variable, - const Number dummy) + const std::vector:: + template solution_curl_type> & + ScratchData::get_curls(const std::string &global_vector_name, + const Extractor & variable, + const Number dummy) { const std::vector &independent_local_dofs = get_local_dof_values(global_vector_name, dummy); @@ -1344,8 +1316,8 @@ namespace MeshWorker // Now build the return type using RetType = - std::vector::curl_type>; + std::vector:: + template solution_curl_type>; RetType &ret = internal_data_storage.template get_or_add_object_with_name( diff --git a/source/fe/fe_values.cc b/source/fe/fe_values.cc index fb39b05f08..e2391ee2c8 100644 --- a/source/fe/fe_values.cc +++ b/source/fe/fe_values.cc @@ -525,17 +525,18 @@ namespace FEValuesViews const ArrayView & dof_values, const Table<2, dealii::Tensor<2, spacedim>> &shape_hessians, const std::vector::ShapeFunctionData> - & shape_function_data, - std::vector::template OutputType< - Number>::laplacian_type> &laplacians) + &shape_function_data, + std::vector:: + template solution_laplacian_type> &laplacians) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = laplacians.size(); - std::fill(laplacians.begin(), - laplacians.end(), - typename Scalar::template OutputType< - Number>::laplacian_type()); + std::fill( + laplacians.begin(), + laplacians.end(), + typename Scalar::template solution_laplacian_type()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -763,17 +764,18 @@ namespace FEValuesViews const ArrayView & dof_values, const Table<2, dealii::Tensor<1, spacedim>> &shape_gradients, const std::vector::ShapeFunctionData> - & shape_function_data, - std::vector::template OutputType< - Number>::divergence_type> &divergences) + &shape_function_data, + std::vector:: + template solution_divergence_type> &divergences) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = divergences.size(); - std::fill(divergences.begin(), - divergences.end(), - typename Vector::template OutputType< - Number>::divergence_type()); + std::fill( + divergences.begin(), + divergences.end(), + typename Vector::template solution_divergence_type()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -1087,17 +1089,18 @@ namespace FEValuesViews const ArrayView & dof_values, const Table<2, dealii::Tensor<2, spacedim>> &shape_hessians, const std::vector::ShapeFunctionData> - & shape_function_data, - std::vector::template OutputType< - Number>::laplacian_type> &laplacians) + &shape_function_data, + std::vector:: + template solution_laplacian_type> &laplacians) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = laplacians.size(); - std::fill(laplacians.begin(), - laplacians.end(), - typename Vector::template OutputType< - Number>::laplacian_type()); + std::fill( + laplacians.begin(), + laplacians.end(), + typename Vector::template solution_laplacian_type()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -1228,15 +1231,15 @@ namespace FEValuesViews typename SymmetricTensor<2, dim, spacedim>::ShapeFunctionData> &shape_function_data, std::vector:: - template OutputType::divergence_type> &divergences) + template solution_divergence_type> &divergences) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = divergences.size(); std::fill(divergences.begin(), divergences.end(), - typename SymmetricTensor<2, dim, spacedim>::template OutputType< - Number>::divergence_type()); + typename SymmetricTensor<2, dim, spacedim>:: + template solution_divergence_type()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -1402,17 +1405,18 @@ namespace FEValuesViews const ArrayView & dof_values, const Table<2, dealii::Tensor<1, spacedim>> &shape_gradients, const std::vector::ShapeFunctionData> - & shape_function_data, - std::vector::template OutputType< - Number>::divergence_type> &divergences) + &shape_function_data, + std::vector:: + template solution_divergence_type> &divergences) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = divergences.size(); - std::fill(divergences.begin(), - divergences.end(), - typename Tensor<2, dim, spacedim>::template OutputType< - Number>::divergence_type()); + std::fill( + divergences.begin(), + divergences.end(), + typename Tensor<2, dim, spacedim>::template solution_divergence_type< + Number>()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -1471,17 +1475,18 @@ namespace FEValuesViews const ArrayView & dof_values, const Table<2, dealii::Tensor<1, spacedim>> &shape_gradients, const std::vector::ShapeFunctionData> - & shape_function_data, - std::vector::template OutputType< - Number>::gradient_type> &gradients) + &shape_function_data, + std::vector:: + template solution_gradient_type> &gradients) { const unsigned int dofs_per_cell = dof_values.size(); const unsigned int n_quadrature_points = gradients.size(); - std::fill(gradients.begin(), - gradients.end(), - typename Tensor<2, dim, spacedim>::template OutputType< - Number>::gradient_type()); + std::fill( + gradients.begin(), + gradients.end(), + typename Tensor<2, dim, spacedim>::template solution_gradient_type< + Number>()); for (unsigned int shape_function = 0; shape_function < dofs_per_cell; ++shape_function) @@ -1541,9 +1546,8 @@ namespace FEValuesViews void Scalar::get_function_values( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1573,9 +1577,8 @@ namespace FEValuesViews void Scalar::get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1598,8 +1601,7 @@ namespace FEValuesViews void Scalar::get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -1629,8 +1631,7 @@ namespace FEValuesViews void Scalar::get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -1654,8 +1655,7 @@ namespace FEValuesViews void Scalar::get_function_hessians( const InputVector &fe_function, - std::vector::type> + std::vector> &hessians) const { Assert(fe_values->update_flags & update_hessians, @@ -1685,8 +1685,7 @@ namespace FEValuesViews void Scalar::get_function_hessians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::hessian_type> + std::vector> &hessians) const { Assert(fe_values->update_flags & update_hessians, @@ -1710,8 +1709,7 @@ namespace FEValuesViews void Scalar::get_function_laplacians( const InputVector &fe_function, - std::vector< - typename ProductType::type> + std::vector> &laplacians) const { Assert(fe_values->update_flags & update_hessians, @@ -1741,8 +1739,7 @@ namespace FEValuesViews void Scalar::get_function_laplacians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::laplacian_type> + std::vector> &laplacians) const { Assert(fe_values->update_flags & update_hessians, @@ -1766,8 +1763,8 @@ namespace FEValuesViews void Scalar::get_function_third_derivatives( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_third_derivative_type> &third_derivatives) const { Assert(fe_values->update_flags & update_3rd_derivatives, @@ -1796,9 +1793,10 @@ namespace FEValuesViews template void Scalar::get_function_third_derivatives_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - third_derivative_type> &third_derivatives) const + const InputVector &dof_values, + std::vector< + solution_third_derivative_type> + &third_derivatives) const { Assert(fe_values->update_flags & update_3rd_derivatives, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1821,9 +1819,8 @@ namespace FEValuesViews void Vector::get_function_values( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1852,9 +1849,8 @@ namespace FEValuesViews void Vector::get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1877,8 +1873,7 @@ namespace FEValuesViews void Vector::get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -1908,8 +1903,7 @@ namespace FEValuesViews void Vector::get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -1933,8 +1927,8 @@ namespace FEValuesViews void Vector::get_function_symmetric_gradients( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_symmetric_gradient_type> &symmetric_gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -1963,9 +1957,10 @@ namespace FEValuesViews template void Vector::get_function_symmetric_gradients_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - symmetric_gradient_type> &symmetric_gradients) const + const InputVector &dof_values, + std::vector< + solution_symmetric_gradient_type> + &symmetric_gradients) const { Assert(fe_values->update_flags & update_gradients, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -1988,8 +1983,7 @@ namespace FEValuesViews void Vector::get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2020,8 +2014,7 @@ namespace FEValuesViews void Vector::get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2045,9 +2038,8 @@ namespace FEValuesViews void Vector::get_function_curls( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &curls) const + std::vector> &curls) + const { Assert(fe_values->update_flags & update_gradients, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2076,8 +2068,7 @@ namespace FEValuesViews void Vector::get_function_curls_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::curl_type> &curls) + std::vector> &curls) const { Assert(fe_values->update_flags & update_gradients, @@ -2101,8 +2092,7 @@ namespace FEValuesViews void Vector::get_function_hessians( const InputVector &fe_function, - std::vector::type> + std::vector> &hessians) const { Assert(fe_values->update_flags & update_hessians, @@ -2132,8 +2122,7 @@ namespace FEValuesViews void Vector::get_function_hessians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::hessian_type> + std::vector> &hessians) const { Assert(fe_values->update_flags & update_hessians, @@ -2157,8 +2146,7 @@ namespace FEValuesViews void Vector::get_function_laplacians( const InputVector &fe_function, - std::vector< - typename ProductType::type> + std::vector> &laplacians) const { Assert(fe_values->update_flags & update_hessians, @@ -2193,8 +2181,7 @@ namespace FEValuesViews void Vector::get_function_laplacians_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::laplacian_type> + std::vector> &laplacians) const { Assert(fe_values->update_flags & update_hessians, @@ -2221,8 +2208,8 @@ namespace FEValuesViews void Vector::get_function_third_derivatives( const InputVector &fe_function, - std::vector::type> + std::vector< + solution_third_derivative_type> &third_derivatives) const { Assert(fe_values->update_flags & update_3rd_derivatives, @@ -2251,9 +2238,10 @@ namespace FEValuesViews template void Vector::get_function_third_derivatives_from_local_dof_values( - const InputVector & dof_values, - std::vector:: - third_derivative_type> &third_derivatives) const + const InputVector &dof_values, + std::vector< + solution_third_derivative_type> + &third_derivatives) const { Assert(fe_values->update_flags & update_3rd_derivatives, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2276,9 +2264,8 @@ namespace FEValuesViews void SymmetricTensor<2, dim, spacedim>::get_function_values( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2307,9 +2294,8 @@ namespace FEValuesViews void SymmetricTensor<2, dim, spacedim>::get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2332,8 +2318,7 @@ namespace FEValuesViews void SymmetricTensor<2, dim, spacedim>::get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2365,8 +2350,7 @@ namespace FEValuesViews SymmetricTensor<2, dim, spacedim>:: get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2390,9 +2374,8 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_values( const InputVector &fe_function, - std::vector< - typename ProductType::type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2421,9 +2404,8 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_values_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::value_type> - &values) const + std::vector> &values) + const { Assert(fe_values->update_flags & update_values, (typename FEValuesBase::ExcAccessToUninitializedField( @@ -2446,8 +2428,7 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_divergences( const InputVector &fe_function, - std::vector::type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2478,8 +2459,7 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_divergences_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::divergence_type> + std::vector> &divergences) const { Assert(fe_values->update_flags & update_gradients, @@ -2503,8 +2483,7 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_gradients( const InputVector &fe_function, - std::vector::type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, @@ -2535,8 +2514,7 @@ namespace FEValuesViews void Tensor<2, dim, spacedim>::get_function_gradients_from_local_dof_values( const InputVector &dof_values, - std::vector< - typename OutputType::gradient_type> + std::vector> &gradients) const { Assert(fe_values->update_flags & update_gradients, diff --git a/source/fe/fe_values.inst.in b/source/fe/fe_values.inst.in index 675d2bbb89..fceba81020 100644 --- a/source/fe/fe_values.inst.in +++ b/source/fe/fe_values.inst.in @@ -249,119 +249,110 @@ for (VEC : GENERAL_CONTAINER_TYPES; Number : ALL_SCALAR_TYPES; template void FEValuesViews::Scalar:: get_function_values_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Scalar:: get_function_gradients_from_local_dof_values>( - const VEC &, - std::vector::gradient_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Scalar:: get_function_hessians_from_local_dof_values>( - const VEC &, - std::vector::hessian_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Scalar:: get_function_laplacians_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Scalar:: get_function_third_derivatives_from_local_dof_values>( const VEC &, - std::vector::third_derivative_type> &) - const; + std::vector> &) const; template void FEValuesViews::Vector:: get_function_values_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Vector:: get_function_gradients_from_local_dof_values>( - const VEC &, - std::vector::gradient_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Vector:: get_function_symmetric_gradients_from_local_dof_values>( const VEC &, - std::vector::symmetric_gradient_type> &) - const; + std::vector> &) const; template void FEValuesViews::Vector:: get_function_divergences_from_local_dof_values>( - const VEC &, - std::vector::divergence_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Vector:: get_function_curls_from_local_dof_values>( - const VEC &, - std::vector::curl_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Vector:: get_function_hessians_from_local_dof_values>( - const VEC &, - std::vector::hessian_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Vector:: get_function_laplacians_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Vector:: get_function_third_derivatives_from_local_dof_values>( const VEC &, - std::vector::third_derivative_type> &) - const; + std::vector> &) const; - template void FEValuesViews:: - SymmetricTensor<2, deal_II_dimension, deal_II_space_dimension>:: - get_function_values_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + template void FEValuesViews::SymmetricTensor<2, + deal_II_dimension, + deal_II_space_dimension>:: + get_function_values_from_local_dof_values>( + const VEC &, std::vector> &) const; template void FEValuesViews:: SymmetricTensor<2, deal_II_dimension, deal_II_space_dimension>:: get_function_divergences_from_local_dof_values>( - const VEC &, - std::vector::divergence_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Tensor<2, deal_II_dimension, deal_II_space_dimension>:: get_function_values_from_local_dof_values>( - const VEC &, - std::vector::value_type> &) const; + const VEC &, std::vector> &) const; template void FEValuesViews::Tensor<2, deal_II_dimension, deal_II_space_dimension>:: get_function_divergences_from_local_dof_values>( - const VEC &, - std::vector::divergence_type> &) const; + const VEC &, std::vector> &) + const; template void FEValuesViews::Tensor<2, deal_II_dimension, deal_II_space_dimension>:: get_function_gradients_from_local_dof_values>( - const VEC &, - std::vector::gradient_type> &) const; + const VEC &, std::vector> &) + const; #endif }