Skip to content

Changed Flow into Turb for limiter Test#403

Closed
jschuele wants to merge 7 commits intosu2code:developfrom
jschuele:js_turb
Closed

Changed Flow into Turb for limiter Test#403
jschuele wants to merge 7 commits intosu2code:developfrom
jschuele:js_turb

Conversation

@jschuele
Copy link

No description provided.

@paulhzhang
Copy link
Contributor

Hi Josef, I agree with you on this modification. The turbulence model can use a different MUSCL scheme coefficient and limiter, as shown in "Grid Convergence for Turbulence Flows (Invited)" by Boris Diskin at el. The second-order extrapolation in turbulence convection can work well without limiting sometimes.

@talbring
Copy link
Member

talbring commented Jul 2, 2017

Nice catch Josef, thanks. I just removed the call to the flow limiter computation altogether.

@economon
Copy link
Member

economon commented Jul 6, 2017

This could be a nice way to save some computational cost... but, the flow limiter may be needed even when the turb. limiter is not in the reconstruction, here:

if (limiter) {
Limiter_i = solver_container[FLOW_SOL]->node[iPoint]->GetLimiter_Primitive();
Limiter_j = solver_container[FLOW_SOL]->node[jPoint]->GetLimiter_Primitive();
}

So, we should double-check in the regressions that the flow limiter is always well computed (during the previous flow iteration) when we perform an iteration of the turb. model.

@talbring
Copy link
Member

But is it really necessary to use the values of the current iteration instead of the values from the previous ? I dont think that it has any benefits.

@economon
Copy link
Member

Agreed, I am not too concerned about reusing the limiter values from a previous iteration. Though it appears to me that we are using the same "limiter" boolean for checking whether we need both the mean flow and turbulence limiters in that routine.. perhaps we can separate that into flow and turbulence limiter booleans so they can operate independently (we have GetSpatialOrder_Flow() and GetSpatialOrder_Turb()).

@talbring
Copy link
Member

I used already GetSpatialOrder_Turb() there, what do you mean exactly ?

@economon
Copy link
Member

@talbring : I was referring to the "limiter" boolean that it use in Upwind_Residual() in solver_direct_turbulent.cpp.

@talbring
Copy link
Member

Ok now I think I understand what you mean. What do you think of my recent change then ?

@economon
Copy link
Member

These changes are still related to the preprocessing.. I was thinking more about the actual reconstructions in the Upwind_Residual() routine. We have booleans that should be brought in line there too.

@talbring
Copy link
Member

talbring commented Aug 1, 2017

Isn't that already handled by the call to SetGlobalParam() in CFluidIteration::Iterate() ?

@ghost
Copy link

ghost commented Aug 6, 2017

Hi,

Unfortunately the logic is quite confuse but I think that part of the code is correct, SU2 just want to use an updated version of the limiter for the turbulence variables

if (config->GetSpatialOrder() == SECOND_ORDER_LIMITER) SetSolution_Limiter(geometry, config);
and for the flow variables

if (limiter_flow) solver_container[FLOW_SOL]->SetPrimitive_Limiter(geometry, config);

Another story is if that is really worth it and we should use the values from the previous iteration. Remember that we compute the turbulence variables after updating the flow variables... so... from the theoretical point of view it makes sense to have an updated value for the limiter. However this is not 100% accurate because we are using old gradients to compute the new limiter. Anyway, the spirit is to have a limiter computed using the latest flow variables value and we should probably check that the gradients have been also updated before computing limiters again for the turbulence equation.

What I don't really like is Upwind_Residual

bool limiter = (config->GetSpatialOrder() == SECOND_ORDER_LIMITER);

in this case we are using only the info about limiters or not in the turbulence model... so if, for example, we are using JST for the mean flow equations then the limiter for the flow equations have never been computed and

        if (limiter) {
          FlowPrimVar_i[iVar] = V_i[iVar] + Limiter_i[iVar]*Project_Grad_i;
          FlowPrimVar_j[iVar] = V_j[iVar] + Limiter_j[iVar]*Project_Grad_j;
        }

Is incorrect because Limiter_i has never been computed (if we are lucky maybe we are using a zero value by default... but we need to check). In other words, I think this require a deeper look. @jschuele, could you please think on this and submit a new pull request (I'll close this one). If you don't have time, please let me know and we really appreciate if you could add this situation as an Issue. We'll try to fix this problem as soon as possible.

By the way the initial confusion was that the SpatialOrder in

config->GetSpatialOrder()
is updated in

void CConfig::SetGlobalParam(unsigned short val_solver,
                             unsigned short val_system,
                             unsigned long val_extiter) {

so it will change during runtime and as we are solving the turbulence model

config->GetSpatialOrder() is the same as config->GetSpatialOrder_Turb()

in
void CTurbSASolver::Preprocessing(CGeometry *geometry, CSolver **solver_container, CConfig *config, unsigned short iMesh, unsigned short iRKStep, unsigned short RunTime_EqSystem, bool Output) {

@jschuele and @paulhzhang thanks a lot for working with us on SU2!
Best,
Francisco

@ghost ghost closed this Aug 6, 2017
@jschuele
Copy link
Author

jschuele commented Aug 7, 2017

I don't get what you asked me to do. I just debugged the turb_oneram case and recognized that limiters were computed twice but only used once within each iteration.
If you give me more hints what to do, I can work on this.

This pull request was closed.
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

4 participants