Skip to content

Commit 81a0620

Browse files
committed
add option specific to turbulence for accurate jacobians
1 parent 1148082 commit 81a0620

4 files changed

Lines changed: 14 additions & 1 deletion

File tree

Common/include/CConfig.hpp

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -590,6 +590,7 @@ class CConfig {
590590
MUSCL_AdjTurb; /*!< \brief MUSCL scheme for the adj turbulence equations.*/
591591
bool MUSCL_Species; /*!< \brief MUSCL scheme for the species equations.*/
592592
bool Use_Accurate_Jacobians; /*!< \brief Use numerically computed Jacobians for AUSM+up(2) and SLAU(2). */
593+
bool Use_Accurate_Turb_Jacobians; /*!< \brief Use numerically computed Jacobians for standard SA turbulence model. */
593594
bool EulerPersson; /*!< \brief Boolean to determine whether this is an Euler simulation with Persson shock capturing. */
594595
bool FSI_Problem = false,/*!< \brief Boolean to determine whether the simulation is FSI or not. */
595596
Multizone_Problem; /*!< \brief Boolean to determine whether we are solving a multizone problem. */
@@ -4499,6 +4500,12 @@ class CConfig {
44994500
*/
45004501
bool GetUse_Accurate_Jacobians(void) const { return Use_Accurate_Jacobians; }
45014502

4503+
/*!
4504+
* \brief Get whether to "Use Accurate Jacobians" for Standard SA turbulence model.
4505+
* \return yes/no.
4506+
*/
4507+
bool GetUse_Accurate_Turb_Jacobians(void) const { return Use_Accurate_Turb_Jacobians; }
4508+
45024509
/*!
45034510
* \brief Get the kind of integration scheme (explicit or implicit)
45044511
* for the flow equations.

Common/src/CConfig.cpp

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1996,6 +1996,8 @@ void CConfig::SetConfig_Options() {
19961996
/*!\brief CONV_NUM_METHOD_TURB
19971997
* \n DESCRIPTION: Convective numerical method \ingroup Config*/
19981998
addConvectOption("CONV_NUM_METHOD_TURB", Kind_ConvNumScheme_Turb, Kind_Centered_Turb, Kind_Upwind_Turb);
1999+
/*!\brief USE_ACCURATE_TURB_JACOBIANS \n DESCRIPTION: Use numerically computed Jacobians for Standard SA model \ingroup Config*/
2000+
addBoolOption("USE_ACCURATE_TURB_JACOBIANS", Use_Accurate_Turb_Jacobians, false);
19992001

20002002
/*!\brief MUSCL_ADJTURB \n DESCRIPTION: Check if the MUSCL scheme should be used \ingroup Config*/
20012003
addBoolOption("MUSCL_ADJTURB", MUSCL_AdjTurb, false);

SU2_CFD/include/numerics/turbulent/turb_diffusion.hpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -93,7 +93,7 @@ class CAvgGrad_TurbSA final : public CAvgGrad_Scalar<FlowIndices> {
9393

9494
/*--- For Jacobians -> Use of TSL approx. to compute derivatives of the gradients ---*/
9595
su2double diffusion_coefficient = term_1 - term_2;
96-
bool UseExactJacobians = config->GetUse_Accurate_Jacobians();
96+
bool UseExactJacobians = config->GetUse_Accurate_Turb_Jacobians();
9797

9898
if (implicit) {
9999
if (UseExactJacobians) {

config_template.cfg

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1849,6 +1849,10 @@ CFL_REDUCTION_TURB= 1.0
18491849
LOWER_LIMIT_K_FACTOR= 1e-15
18501850
LOWER_LIMIT_OMEGA_FACTOR= 1e-5
18511851
1852+
% Use numerically computed exact Jacobians for standard SA turbulence model
1853+
% Slower per iteration but potentialy more stable and capable of higher CFL
1854+
USE_ACCURATE_FLUX_JACOBIANS= NO
1855+
%
18521856
% --------------------- HEAT NUMERICAL METHOD DEFINITION ----------------------%
18531857
%
18541858
% Value of the thermal diffusivity

0 commit comments

Comments
 (0)