Skip to content
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

Hybrid Parallel AD (Part 2/?) #1284

Merged
merged 31 commits into from
May 29, 2021
Merged
Show file tree
Hide file tree
Changes from 1 commit
Commits
Show all changes
31 commits
Select commit Hold shift + click to select a range
f7f1b96
local indices to allow parallel registration
pcarruscag May 9, 2021
20cbe7a
reset inputs as they are extracted
pcarruscag May 9, 2021
e8522c6
fix #1285
pcarruscag May 10, 2021
7f9f647
address #1273
pcarruscag May 10, 2021
2144acd
try to fix restarts
pcarruscag May 13, 2021
7efd040
Merge remote-tracking branch 'upstream/develop' into hybrid_parallel_ad3
pcarruscag May 13, 2021
c4e3eb4
Merge remote-tracking branch 'upstream/feature_dynFSIAD' into hybrid_…
pcarruscag May 13, 2021
0d55909
update regressions
pcarruscag May 17, 2021
0f06d7a
Merge remote-tracking branch 'upstream/add_unstchtcase' into hybrid_p…
pcarruscag May 17, 2021
91503ea
Merge branch 'feature_dynFSIAD' into hybrid_parallel_ad3
pcarruscag May 18, 2021
8c6233a
Move \!Crossterm bool in front of loop for unsteady adjoint extraction.
TobiKattmann May 18, 2021
7422c00
dont extract some terms during cross term evaluation
pcarruscag May 18, 2021
d598915
Merge branch 'hybrid_parallel_ad3' of /~https://github.com/su2code/SU2 …
pcarruscag May 18, 2021
d949cbe
apply BGS relaxation also to velocities
pcarruscag May 18, 2021
b7f5b05
cleanup overzealous SU2_OMP_MASTER
pcarruscag May 19, 2021
dcc5038
Merge remote-tracking branch 'origin/develop' into hybrid_parallel_ad3
TobiKattmann May 24, 2021
05d1ce7
Merge remote-tracking branch 'upstream/develop' into hybrid_parallel_ad3
pcarruscag May 24, 2021
6364cbf
less duplication, better unsteady FSI restarts
pcarruscag May 24, 2021
d2fbc13
fix for turbo
pcarruscag May 24, 2021
cd9b974
Merge branch 'hybrid_parallel_ad3' of /~https://github.com/su2code/SU2 …
TobiKattmann May 24, 2021
319e51e
fix for moving grid disc adj
pcarruscag May 24, 2021
6094cb6
update tests
pcarruscag May 24, 2021
d5715c4
local indices only, no more hidden global "adjointPosition"
pcarruscag May 25, 2021
567064a
fixes
pcarruscag May 25, 2021
ef4f233
proper-er way to reset input
pcarruscag May 26, 2021
9df5d02
small things
pcarruscag May 26, 2021
9a6346f
Merge branch 'hybrid_parallel_ad3' of /~https://github.com/su2code/SU2 …
TobiKattmann May 27, 2021
19d3329
Remove errors from #1281. No auto for su2double with rhs math operati…
TobiKattmann May 27, 2021
e8cee0b
CoDiPack update.
jblueh May 27, 2021
be27e1d
Testing RealReverseIndex.
jblueh May 27, 2021
4593083
Merge branch 'develop' into hybrid_parallel_ad3
pcarruscag May 29, 2021
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
55 changes: 16 additions & 39 deletions SU2_CFD/src/solvers/CDiscAdjFEASolver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -236,64 +236,41 @@ void CDiscAdjFEASolver::RegisterOutput(CGeometry *geometry, CConfig *config){

}

void CDiscAdjFEASolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *config, bool CrossTerm){

const bool dynamic = config->GetTime_Domain();
const bool multizone = config->GetMultizone_Problem();

unsigned short iVar;
unsigned long iPoint;
su2double residual;

su2double Solution[MAXNVAR] = {0.0};

/*--- Set Residuals to zero ---*/

SetResToZero();
void CDiscAdjFEASolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *config, bool CrossTerm) {

/*--- Set the old solution, for multi-zone problems this is done after computing the
* residuals, otherwise the per-zone-residuals do not make sense, as on entry Solution
* contains contributions from other zones but on extraction it does not. ---*/

if(!multizone) nodes->Set_OldSolution();

for (iPoint = 0; iPoint < nPoint; iPoint++){
if (!config->GetMultizone_Problem()) nodes->Set_OldSolution();

/*--- Extract the adjoint solution ---*/
/*--- Extract and store the adjoint solution ---*/

for (auto iPoint = 0ul; iPoint < nPoint; iPoint++) {
su2double Solution[MAXNVAR] = {0.0};
direct_solver->GetNodes()->GetAdjointSolution(iPoint,Solution);

/*--- Store the adjoint solution ---*/

nodes->SetSolution(iPoint,Solution);

}

/*--- Solution for acceleration (u'') and velocity (u') at time n ---*/
if (CrossTerm) return;

if (dynamic){

/*--- NOW: The solution at time n ---*/
for (iPoint = 0; iPoint < nPoint; iPoint++){

/*--- Extract the adjoint solution at time n ---*/
/*--- Extract and store the adjoint solution at time n (including accel. and velocity) ---*/

if (config->GetTime_Domain()) {
for (auto iPoint = 0ul; iPoint < nPoint; iPoint++) {
su2double Solution[MAXNVAR] = {0.0};
direct_solver->GetNodes()->GetAdjointSolution_time_n(iPoint,Solution);

/*--- Store the adjoint solution at time n ---*/

if (!CrossTerm) nodes->Set_Solution_time_n(iPoint,Solution);
nodes->Set_Solution_time_n(iPoint,Solution);
}

}

/*--- TODO: Need to set the MPI solution in the previous TS ---*/

/*--- Set the residuals ---*/

for (iPoint = 0; iPoint < nPointDomain; iPoint++){
for (iVar = 0; iVar < nVar; iVar++){
residual = nodes->GetSolution(iPoint, iVar) - nodes->GetSolution_Old(iPoint, iVar);
SetResToZero();

for (auto iPoint = 0ul; iPoint < nPointDomain; iPoint++) {
for (auto iVar = 0u; iVar < nVar; iVar++){
su2double residual = nodes->GetSolution(iPoint, iVar) - nodes->GetSolution_Old(iPoint, iVar);

Residual_RMS[iVar] += residual*residual;
AddRes_Max(iVar,fabs(residual),geometry->nodes->GetGlobalIndex(iPoint),geometry->nodes->GetCoord(iPoint));
Expand Down
26 changes: 13 additions & 13 deletions SU2_CFD/src/solvers/CDiscAdjSolver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -297,32 +297,29 @@ void CDiscAdjSolver::RegisterOutput(CGeometry *geometry, CConfig *config) {
direct_solver->GetNodes()->RegisterSolution(false);
}

void CDiscAdjSolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *config, bool CrossTerm){
void CDiscAdjSolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *config, bool CrossTerm) {

const bool time_n1_needed = config->GetTime_Marching() == TIME_MARCHING::DT_STEPPING_2ND;
const bool time_n_needed = (config->GetTime_Marching() == TIME_MARCHING::DT_STEPPING_1ST) || time_n1_needed;

const su2double relax = (config->GetInnerIter()==0) ? 1.0 : config->GetRelaxation_Factor_Adjoint();

su2double Solution[MAXNVAR] = {0.0};

/*--- Set Residuals to zero ---*/

SetResToZero();
/*--- Thread-local residual variables. ---*/

su2double resMax[MAXNVAR] = {0.0}, resRMS[MAXNVAR] = {0.0};
const su2double* coordMax[MAXNVAR] = {nullptr};
unsigned long idxMax[MAXNVAR] = {0};

/*--- Set the old solution and compute residuals. ---*/

if(!config->GetMultizone_Problem()) nodes->Set_OldSolution();
if (!config->GetMultizone_Problem()) nodes->Set_OldSolution();

SU2_OMP_FOR_STAT(omp_chunk_size)
for (auto iPoint = 0ul; iPoint < nPoint; iPoint++) {

/*--- Extract the adjoint solution ---*/

su2double Solution[MAXNVAR] = {0.0};
pcarruscag marked this conversation as resolved.
Show resolved Hide resolved
direct_solver->GetNodes()->GetAdjointSolution(iPoint,Solution);
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Nice to have the input-index & output-index storing-strategies consolidated between single and multizone adjoint 👍


/*--- Relax and store the adjoint solution, compute the residuals. ---*/
Expand All @@ -344,6 +341,11 @@ void CDiscAdjSolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *confi
}
END_SU2_OMP_FOR

/*--- Residuals and time_n terms are not needed when evaluating multizone cross terms. ---*/
if (CrossTerm) return;
pcarruscag marked this conversation as resolved.
Show resolved Hide resolved

SetResToZero();
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Could we have an unwanted addition to the residual because of moving the SetResToZero below the CrossTerm-return? The residual should compare the (full) adjoint solution so aren't we missing the External (crossterms + dualTime terms) here anyway. I guess just the dualTime Terms are constant so we can leave them out. Or is it correct now as CrossTerms are extracted first. I have to look that up

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

No, the residual should measure the norm of the residual of the (linear) adjoint equations, otherwise it would measure stagnation instead of convergence.
It just so happens that differences between consecutive fixed point iterations give you that residual (because the fixed point is the same as a Richardson iteration for the right-preconditioned adjoint equations).
The single zone uses differences of the total terms, the multizone uses differences of the product between total terms and the Jacobian of the iteration (it is easier to do it that way), this is equivalent since the right hand side is constant during inner iterations.


/*--- Reduce residual information over all threads in this rank. ---*/
SU2_OMP_CRITICAL
for (auto iVar = 0u; iVar < nVar; iVar++) {
Expand All @@ -365,10 +367,9 @@ void CDiscAdjSolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *confi
if (time_n_needed) {
SU2_OMP_FOR_STAT(omp_chunk_size)
for (auto iPoint = 0ul; iPoint < nPoint; iPoint++) {

su2double Solution[MAXNVAR] = {0.0};
direct_solver->GetNodes()->GetAdjointSolution_time_n(iPoint,Solution);

if (!CrossTerm) nodes->Set_Solution_time_n(iPoint,Solution);
nodes->Set_Solution_time_n(iPoint,Solution);
}
END_SU2_OMP_FOR
}
Expand All @@ -377,10 +378,9 @@ void CDiscAdjSolver::ExtractAdjoint_Solution(CGeometry *geometry, CConfig *confi
if (time_n1_needed) {
SU2_OMP_FOR_STAT(omp_chunk_size)
for (auto iPoint = 0ul; iPoint < nPoint; iPoint++) {

su2double Solution[MAXNVAR] = {0.0};
direct_solver->GetNodes()->GetAdjointSolution_time_n1(iPoint,Solution);

if (!CrossTerm) nodes->Set_Solution_time_n1(iPoint,Solution);
nodes->Set_Solution_time_n1(iPoint,Solution);
}
END_SU2_OMP_FOR
}
Expand Down