@@ -62,7 +62,6 @@ void CIncNSSolver::Preprocessing(CGeometry *geometry, CSolver **solver_container
6262 const bool center = (config->GetKind_ConvNumScheme_Flow () == SPACE_CENTERED);
6363 const bool limiter = (config->GetKind_SlopeLimit_Flow () != NO_LIMITER) && (InnerIter <= config->GetLimiterIter ());
6464 const bool van_albada = (config->GetKind_SlopeLimit_Flow () == VAN_ALBADA_EDGE);
65- const bool energy = config->GetEnergy_Equation ();
6665
6766 /* --- Common preprocessing steps (implemented by CEulerSolver) ---*/
6867
@@ -99,42 +98,49 @@ void CIncNSSolver::Preprocessing(CGeometry *geometry, CSolver **solver_container
9998 ComputeVorticityAndStrainMag<1 >(*config, iMesh);
10099
101100 /* --- Compute recovered pressure and temperature for streamwise periodic flow ---*/
102- if (config->GetKind_Streamwise_Periodic () != NONE) {
101+ if (config->GetKind_Streamwise_Periodic () != NONE)
102+ Compute_Streamwise_Periodic_Recovered_Values (config, geometry, iMesh);
103+ }
103104
104- /* --- Define and initialize helping variables --- */
105- su2double dot_product, Pressure_Recovered, Temperature_Recovered;
105+ void CIncNSSolver::Compute_Streamwise_Periodic_Recovered_Values (CConfig *config, const CGeometry *geometry,
106+ const unsigned short iMesh) {
106107
107- const su2double delta_p = config->GetStreamwise_Periodic_PressureDrop () / config->GetPressure_Ref ();
108+ const bool energy = (config->GetEnergy_Equation () && config->GetStreamwise_Periodic_Temperature ());
109+ const auto InnerIter = config->GetInnerIter ();
108110
109- /* --- Reference node on inlet periodic marker to compute relative distance along periodic translation vector. ---*/
110- const su2double* ReferenceNode = geometry->GetStreamwise_Periodic_RefNode ();
111+ const su2double delta_p = config->GetStreamwise_Periodic_PressureDrop () / config->GetPressure_Ref ();
111112
112- /* --- Compute square of the distance between the 2 periodic surfaces . ---*/
113- const su2double norm2_translation = GeometryToolbox::SquaredNorm (nDim, config-> GetPeriodic_Translation ( 0 ) );
113+ /* --- Reference node on inlet periodic marker to compute relative distance along periodic translation vector . ---*/
114+ const su2double* ReferenceNode = geometry-> GetStreamwise_Periodic_RefNode ( );
114115
115- /* --- Compute recoverd pressure and temperature for all points ---*/
116- for ( auto iPoint = 0ul ; iPoint < nPoint; iPoint++) {
116+ /* --- Compute square of the distance between the 2 periodic surfaces. ---*/
117+ const su2double norm2_translation = GeometryToolbox::SquaredNorm (nDim, config-> GetPeriodic_Translation ( 0 ));
117118
118- /* --- First, compute helping terms based on relative distance (0,l) between periodic markers ---*/
119- dot_product = 0.0 ;
120- for (unsigned short iDim = 0 ; iDim < nDim; iDim++)
121- dot_product += fabs ( (geometry->nodes ->GetCoord (iPoint,iDim) - ReferenceNode[iDim]) * config->GetPeriodic_Translation (0 )[iDim]);
119+ /* --- Compute recoverd pressure and temperature for all points ---*/
120+ SU2_OMP_FOR_STAT (omp_chunk_size)
121+ for (auto iPoint = 0ul ; iPoint < nPoint; iPoint++) {
122122
123- /* --- Second, substract/add correction from reduced pressure/temperature to get recoverd pressure/temperature ---*/
124- Pressure_Recovered = nodes->GetSolution (iPoint, 0 ) - delta_p / norm2_translation * dot_product;
125- nodes->SetStreamwise_Periodic_RecoveredPressure (iPoint, Pressure_Recovered);
123+ /* --- First, compute helping terms based on relative distance (0,l) between periodic markers ---*/
124+ su2double dot_product = 0.0 ;
125+ for (unsigned short iDim = 0 ; iDim < nDim; iDim++)
126+ dot_product += fabs ( (geometry->nodes ->GetCoord (iPoint,iDim) - ReferenceNode[iDim]) * config->GetPeriodic_Translation (0 )[iDim]);
126127
127- /* --- InnerIter > 0 as otherwise MassFlow in the denominator would be zero ---*/
128- if (energy && InnerIter > 0 ) {
129- Temperature_Recovered = nodes->GetSolution (iPoint, nDim+1 );
130- Temperature_Recovered += Streamwise_Periodic_IntegratedHeatFlow / (Streamwise_Periodic_MassFlow * nodes->GetSpecificHeatCp (iPoint) * norm2_translation) * dot_product;
131- nodes->SetStreamwise_Periodic_RecoveredTemperature (iPoint, Temperature_Recovered);
132- }
128+ /* --- Second, substract/add correction from reduced pressure/temperature to get recoverd pressure/temperature ---*/
129+ const su2double Pressure_Recovered = nodes->GetPressure (iPoint) - delta_p / norm2_translation * dot_product;
130+ nodes->SetStreamwise_Periodic_RecoveredPressure (iPoint, Pressure_Recovered);
131+
132+ /* --- InnerIter > 0 as otherwise MassFlow in the denominator would be zero ---*/
133+ if (energy && InnerIter > 0 ) {
134+ su2double Temperature_Recovered = nodes->GetTemperature (iPoint);
135+ Temperature_Recovered += Streamwise_Periodic_IntegratedHeatFlow / (Streamwise_Periodic_MassFlow * nodes->GetSpecificHeatCp (iPoint) * norm2_translation) * dot_product;
136+ nodes->SetStreamwise_Periodic_RecoveredTemperature (iPoint, Temperature_Recovered);
133137 }
138+ } // for iPoint
134139
135- /* --- Compute the integrated Heatflux Q into the domain, and massflow over periodic markers ---*/
136- GetStreamwise_Periodic_Properties (geometry, config, iMesh);
137- } // if streamwise periodic
140+ /* --- Compute the integrated Heatflux Q into the domain, and massflow over periodic markers ---*/
141+ SU2_OMP_MASTER
142+ GetStreamwise_Periodic_Properties (geometry, config, iMesh);
143+ SU2_OMP_BARRIER
138144}
139145
140146void CIncNSSolver::Viscous_Residual (unsigned long iEdge, CGeometry *geometry, CSolver **solver_container,
@@ -265,20 +271,20 @@ void CIncNSSolver::BC_Wall_Generic(const CGeometry *geometry, const CConfig *con
265271 LinSysRes (iPoint, nDim+1 ) -= Wall_HeatFlux*Area;
266272
267273 /* --- With streamwise periodic flow and heatflux walls an additional term is introduced in the boundary formulation ---*/
268- if (streamwise_periodic && streamwise_periodic_temperature) {
274+ if (streamwise_periodic && streamwise_periodic_temperature) {
269275
270- Cp = nodes->GetSpecificHeatCp (iPoint);
271- thermal_conductivity = nodes->GetThermalConductivity (iPoint);
276+ Cp = nodes->GetSpecificHeatCp (iPoint);
277+ thermal_conductivity = nodes->GetThermalConductivity (iPoint);
272278
273- /* --- Scalar factor of the residual contribution ---*/
274- const su2double norm2_translation = GeometryToolbox::SquaredNorm (nDim, config->GetPeriodic_Translation (0 ));
275- scalar_factor = Streamwise_Periodic_IntegratedHeatFlow*thermal_conductivity / (Streamwise_Periodic_MassFlow * Cp * norm2_translation);
279+ /* --- Scalar factor of the residual contribution ---*/
280+ const su2double norm2_translation = GeometryToolbox::SquaredNorm (nDim, config->GetPeriodic_Translation (0 ));
281+ scalar_factor = Streamwise_Periodic_IntegratedHeatFlow*thermal_conductivity / (Streamwise_Periodic_MassFlow * Cp * norm2_translation);
276282
277- /* --- Dot product ---*/
278- dot_product = GeometryToolbox::DotProduct (nDim, config->GetPeriodic_Translation (0 ), Normal);
283+ /* --- Dot product ---*/
284+ dot_product = GeometryToolbox::DotProduct (nDim, config->GetPeriodic_Translation (0 ), Normal);
279285
280- LinSysRes (iPoint, nDim+1 ) += scalar_factor*dot_product;
281- } // if streamwise_periodic
286+ LinSysRes (iPoint, nDim+1 ) += scalar_factor*dot_product;
287+ } // if streamwise_periodic
282288 }
283289 else { // ISOTHERMAL
284290
0 commit comments