Skip to content

Commit fc0f151

Browse files
Fix turbulence restart loading regression by adding legacy index fallback to name-based lookup.
1 parent 1a7b94b commit fc0f151

2 files changed

Lines changed: 48 additions & 1 deletion

File tree

SU2_CFD/src/solvers/CTransLMSolver.cpp

Lines changed: 30 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -515,10 +515,39 @@ void CTransLMSolver::LoadRestart(CGeometry** geometry, CSolver*** solver, CConfi
515515
Read_SU2_Restart_ASCII(geometry[MESH_0], config, restart_filename);
516516
}
517517

518-
/*--- Identify indices for LM transition variables ---*/
518+
/*--- Find indices for LM transition variables ---*/
519519
vector<string> target_fields = {"LM_gamma", "LM_Re_t", "LM_gamma_sep", "LM_gamma_eff"};
520520
vector<int> field_indices = FindFieldIndices(target_fields);
521521

522+
/*--- Fallback: Use legacy offset if name not found ---*/
523+
/*--- nDim + 2 covers Density + Momentum + Energy (Compressible)
524+
and Pressure + Velocity + Temperature (Incompressible) ---*/
525+
const unsigned short flow_offset = geometry[MESH_0]->GetnDim() + 2;
526+
527+
int turb_offset = flow_offset;
528+
switch(config->GetKind_Turb_Model()) {
529+
case TURB_MODEL::SA:
530+
turb_offset += 1;
531+
break;
532+
case TURB_MODEL::SST:
533+
turb_offset += 2;
534+
break;
535+
default: break;
536+
}
537+
538+
for (size_t i = 0; i < target_fields.size(); ++i) {
539+
if (field_indices[i] == -1) {
540+
const int fallback_idx = turb_offset + i;
541+
if (Restart_Vars.size() > 1 && fallback_idx < Restart_Vars[1]) {
542+
field_indices[i] = fallback_idx;
543+
if (rank == MASTER_NODE) {
544+
cout << "WARNING: " << target_fields[i] << " field not found in restart file. "
545+
<< "Using fallback index " << fallback_idx << ".\n";
546+
}
547+
}
548+
}
549+
}
550+
522551
/*--- Warn if any fields are missing (Master node only) ---*/
523552
if (rank == MASTER_NODE) {
524553
for (size_t i = 0; i < target_fields.size(); ++i) {

SU2_CFD/src/solvers/CTurbSolver.cpp

Lines changed: 18 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -134,6 +134,24 @@ void CTurbSolver::LoadRestart(CGeometry** geometry, CSolver*** solver, CConfig*
134134
/*--- Find indices for all turbulence variables at once ---*/
135135
vector<int> field_indices = FindFieldIndices(varNames);
136136

137+
/*--- Fallback: Use legacy offset if name not found ---*/
138+
/*--- nDim + 2 covers Density + Momentum + Energy (Compressible)
139+
and Pressure + Velocity + Temperature (Incompressible) ---*/
140+
const unsigned short flow_offset = geometry[MESH_0]->GetnDim() + 2;
141+
142+
for (unsigned short iVar = 0; iVar < nVar; ++iVar) {
143+
if (field_indices[iVar] == -1 && !varNames[iVar].empty()) {
144+
const int fallback_idx = flow_offset + iVar;
145+
if (Restart_Vars.size() > 1 && fallback_idx < Restart_Vars[1]) {
146+
field_indices[iVar] = fallback_idx;
147+
if (rank == MASTER_NODE) {
148+
cout << "WARNING: " << varNames[iVar] << " field not found in restart file. "
149+
<< "Using fallback index " << fallback_idx << ".\n";
150+
}
151+
}
152+
}
153+
}
154+
137155
/*--- Warn if any fields are missing (Master node only) ---*/
138156
if (rank == MASTER_NODE) {
139157
for (unsigned short iVar = 0; iVar < nVar; ++iVar) {

0 commit comments

Comments
 (0)