-
Notifications
You must be signed in to change notification settings - Fork 708
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Work on ShapeInfo::is_supported() #14953
Merged
Merged
Changes from all commits
Commits
File filter
Filter by extension
Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -30,6 +30,7 @@ | |
#include <deal.II/fe/fe.h> | ||
#include <deal.II/fe/fe_dgp.h> | ||
#include <deal.II/fe/fe_dgq.h> | ||
#include <deal.II/fe/fe_nothing.h> | ||
#include <deal.II/fe/fe_poly.h> | ||
#include <deal.II/fe/fe_pyramid_p.h> | ||
#include <deal.II/fe/fe_q.h> | ||
|
@@ -397,14 +398,12 @@ namespace internal | |
return; | ||
} | ||
else if (quad_in.is_tensor_product() == false || | ||
dynamic_cast<const FE_SimplexP<dim> *>( | ||
&fe_in.base_element(base_element_number)) || | ||
dynamic_cast<const FE_SimplexDGP<dim> *>( | ||
&fe_in.base_element(base_element_number)) || | ||
dynamic_cast<const FE_WedgeP<dim> *>( | ||
&fe_in.base_element(base_element_number)) || | ||
dynamic_cast<const FE_PyramidP<dim> *>( | ||
&fe_in.base_element(base_element_number))) | ||
dynamic_cast<const FE_SimplexPoly<dim, dim> *>( | ||
&fe_in.base_element(base_element_number)) != nullptr || | ||
dynamic_cast<const FE_WedgePoly<dim, dim> *>( | ||
&fe_in.base_element(base_element_number)) != nullptr || | ||
dynamic_cast<const FE_PyramidPoly<dim, dim> *>( | ||
&fe_in.base_element(base_element_number)) != nullptr) | ||
Comment on lines
+401
to
+406
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. This was not correct and not an issue, since the first condition normally is true. |
||
{ | ||
// specialization for arbitrary finite elements and quadrature rules | ||
// as needed in the context, e.g., of simplices | ||
|
@@ -1238,6 +1237,9 @@ namespace internal | |
if (dim != spacedim) | ||
return false; | ||
|
||
if (dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe)) | ||
return true; | ||
|
||
for (unsigned int base = 0; base < fe.n_base_elements(); ++base) | ||
{ | ||
const FiniteElement<dim, spacedim> *fe_ptr = &(fe.base_element(base)); | ||
|
@@ -1251,13 +1253,17 @@ namespace internal | |
dynamic_cast<const FE_Poly<dim, spacedim> *>(fe_ptr); | ||
// Simplices are a special case since the polynomial family is not | ||
// indicative of their support | ||
if (dynamic_cast<const FE_SimplexP<dim> *>(fe_poly_ptr) || | ||
dynamic_cast<const FE_SimplexDGP<dim> *>(fe_poly_ptr) || | ||
dynamic_cast<const FE_WedgeP<dim> *>(fe_poly_ptr) || | ||
dynamic_cast<const FE_PyramidP<dim> *>(fe_poly_ptr)) | ||
if (dynamic_cast<const FE_SimplexPoly<dim, dim> *>(fe_poly_ptr) != | ||
nullptr || | ||
dynamic_cast<const FE_WedgePoly<dim, dim> *>(fe_poly_ptr) != | ||
nullptr || | ||
dynamic_cast<const FE_PyramidPoly<dim, dim> *>(fe_poly_ptr) != | ||
nullptr) | ||
return true; | ||
|
||
if (dynamic_cast<const TensorProductPolynomials<dim> *>( | ||
if (dynamic_cast<const TensorProductPolynomials< | ||
dim, | ||
Polynomials::Polynomial<double>> *>( | ||
&fe_poly_ptr->get_poly_space()) == nullptr && | ||
dynamic_cast<const TensorProductPolynomials< | ||
dim, | ||
|
@@ -1269,6 +1275,9 @@ namespace internal | |
nullptr) | ||
return false; | ||
} | ||
else if (dynamic_cast<const FE_Nothing<dim, spacedim> *>(fe_ptr) != | ||
nullptr) | ||
return true; | ||
else | ||
return false; | ||
} | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -16,7 +16,7 @@ DEAL::FE_FaceP<2>(1) supported by MatrixFree: false | |
DEAL::FE_FaceQ<2>(1) supported by MatrixFree: false | ||
DEAL::FE_P1NC supported by MatrixFree: false | ||
DEAL::FE_Nedelec<2>(1) supported by MatrixFree: false | ||
DEAL::FE_Nothing<2>() supported by MatrixFree: false | ||
DEAL::FE_Nothing<2>() supported by MatrixFree: true | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. 2 new supported elements :) |
||
DEAL::FE_Q<1,2>(1) supported by MatrixFree: false | ||
DEAL::FE_Q<2>(1) supported by MatrixFree: true | ||
DEAL::FE_Q_Bubbles<2>(1) supported by MatrixFree: false | ||
|
@@ -25,5 +25,5 @@ DEAL::FE_Q_Hierarchical<2>(1) supported by MatrixFree: true | |
DEAL::FE_Q_iso_Q1<2>(1) supported by MatrixFree: true | ||
DEAL::FE_RannacherTurek<2>(0, 2) supported by MatrixFree: false | ||
DEAL::FE_RaviartThomas<2>(1) supported by MatrixFree: false | ||
DEAL::FE_RaviartThomasNodal<2>(1) supported by MatrixFree: false | ||
DEAL::FE_RaviartThomasNodal<2>(1) supported by MatrixFree: true | ||
DEAL::FE_TraceQ<2>(1) supported by MatrixFree: false |
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
@kronbichler Have I missed any finite elements?