-
Notifications
You must be signed in to change notification settings - Fork 15.3k
[flang][OpenMP] Tolerate compiler directives in loop constructs #169346
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
Changes from all commits
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -796,6 +796,28 @@ static void processTileSizesFromOpenMPConstruct( | |
| } | ||
| } | ||
|
|
||
| static pft::Evaluation *getNestedDoConstruct(pft::Evaluation &eval) { | ||
| for (pft::Evaluation &nested : eval.getNestedEvaluations()) { | ||
| // In an OpenMPConstruct there can be compiler directives: | ||
| // 1 <<OpenMPConstruct>> | ||
| // 2 CompilerDirective: !unroll | ||
| // <<DoConstruct>> -> 8 | ||
| if (nested.getIf<parser::CompilerDirective>()) | ||
| continue; | ||
| // Within a DoConstruct, there can be compiler directives, plus | ||
| // there is a DoStmt before the body: | ||
| // <<DoConstruct>> -> 8 | ||
| // 3 NonLabelDoStmt -> 7: do i = 1, n | ||
| // <<DoConstruct>> -> 7 | ||
| if (nested.getIf<parser::NonLabelDoStmt>()) | ||
| continue; | ||
| assert(nested.getIf<parser::DoConstruct>() && | ||
| "Unexpected construct in the nested evaluations"); | ||
| return &nested; | ||
| } | ||
| llvm_unreachable("Expected do loop to be in the nested evaluations"); | ||
| } | ||
|
Comment on lines
+799
to
+819
Contributor
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. Will this work with loop generating OpenMP constructs?
Contributor
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. If
Contributor
Author
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. We crash on this testcase (we crashed on it before my loop nest parser work). The lowering code does not support that at the moment.
Contributor
Author
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.
|
||
|
|
||
| /// Populates the sizes vector with values if the given OpenMPConstruct | ||
| /// contains a loop construct with an inner tiling construct. | ||
| void collectTileSizesFromOpenMPConstruct( | ||
|
|
@@ -818,7 +840,7 @@ int64_t collectLoopRelatedInfo( | |
| int64_t numCollapse = 1; | ||
|
|
||
| // Collect the loops to collapse. | ||
| lower::pft::Evaluation *doConstructEval = &eval.getFirstNestedEvaluation(); | ||
| lower::pft::Evaluation *doConstructEval = getNestedDoConstruct(eval); | ||
| if (doConstructEval->getIf<parser::DoConstruct>()->IsDoConcurrent()) { | ||
| TODO(currentLocation, "Do Concurrent in Worksharing loop construct"); | ||
| } | ||
|
|
@@ -844,7 +866,7 @@ void collectLoopRelatedInfo( | |
| fir::FirOpBuilder &firOpBuilder = converter.getFirOpBuilder(); | ||
|
|
||
| // Collect the loops to collapse. | ||
| lower::pft::Evaluation *doConstructEval = &eval.getFirstNestedEvaluation(); | ||
| lower::pft::Evaluation *doConstructEval = getNestedDoConstruct(eval); | ||
| if (doConstructEval->getIf<parser::DoConstruct>()->IsDoConcurrent()) { | ||
| TODO(currentLocation, "Do Concurrent in Worksharing loop construct"); | ||
| } | ||
|
|
@@ -885,9 +907,8 @@ void collectLoopRelatedInfo( | |
| iv.push_back(bounds->name.thing.symbol); | ||
| loopVarTypeSize = std::max(loopVarTypeSize, | ||
| bounds->name.thing.symbol->GetUltimate().size()); | ||
| collapseValue--; | ||
| doConstructEval = | ||
| &*std::next(doConstructEval->getNestedEvaluations().begin()); | ||
| if (--collapseValue) | ||
| doConstructEval = getNestedDoConstruct(*doConstructEval); | ||
| } while (collapseValue > 0); | ||
|
|
||
| convertLoopBounds(converter, currentLocation, result, loopVarTypeSize); | ||
|
|
||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,21 @@ | ||
| !RUN: %flang_fc1 -emit-hlfir -fopenmp -fopenmp-version=60 %s -o - | FileCheck %s | ||
|
|
||
| ! Check that this compiles successfully, but not rely on any specific output. | ||
|
|
||
| !CHECK: omp.parallel | ||
|
|
||
| program omp_cdir_crash | ||
| implicit none | ||
| integer, parameter :: n = 10 | ||
| real :: a(n) | ||
| integer :: i | ||
|
|
||
| !$omp parallel do | ||
| !dir$ unroll | ||
| do i = 1, n | ||
| a(i) = real(i) | ||
| end do | ||
| !$omp end parallel do | ||
|
|
||
| print *, 'a(1)=', a(1), ' a(n)=', a(n) | ||
| end program omp_cdir_crash |
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.
Is there a typo here?
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.
I can't see it... 😐 Is there something strange about the wording?
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.
Ok, now I see an improper use of a comma.
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.
Sorry Kryztof, I was wondering whether you were mentioning the compiler directives point twice (once above as part of the skipping the compiler directive and also here).
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.
Yes, the directives could appear in both constructs. I wanted to mention it because the prior code assumed that DoStmt would be the second nested evaluation.