-
Notifications
You must be signed in to change notification settings - Fork 61
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
Tracer flow rate at boundaries #1197
Conversation
d560eee
to
f7c5794
Compare
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Good work! I have a few comments :)
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Really nice! Nothing much to say, only clarifications for my understanding.
@@ -1655,24 +1655,11 @@ CFDDEMSolver<dim>::solve() | |||
while (this->simulation_control->integrate()) | |||
{ | |||
this->simulation_control->print_progression(this->pcout); | |||
bool refinement_step; | |||
if (this->simulation_parameters.mesh_adaptation.refinement_at_frequency) |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Just curious to know why is those if conditions are removed from the solvers and fem-dem?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
It's because in the update_boundary_conditions
functions, the only condition to respect was whether boundary conditions are time dependent. This is already very restrictive and not expensive to check, so to simplify the code I remove the refinement_step
boolean here. It had no effect on any test. I feel like the previous condition was difficult to understand and didn't bring much.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
If it does the same I like it! It is indeed much clearer :)
4d3d401
to
4c9b8da
Compare
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Nice work! Sorry for the delay in my review!
I think your application test is enough, but a unit test could be also nice :)
@@ -5,6 +5,12 @@ The format is based on [Keep a Changelog](http://keepachangelog.com/). | |||
|
|||
## [Master] - 2024-08-05 | |||
|
|||
### Added | |||
|
|||
- MINOR A post-processing feature for the tracer flow rate through boundaries is added. This allows to produce data useful for residence time distribution analyses. [#1197](https://github.com/chaos-polymtl/lethe/pull/1197) |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I think you should add the fix on the bc update
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
one small correction left
Nice work! :)
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
One small refactoring left and this should be good to go.
source/solvers/tracer.cc
Outdated
// Fill table | ||
this->tracer_flow_rate_table.add_value( | ||
"time", this->simulation_control->get_current_time()); | ||
this->tracer_flow_rate_table.set_precision("time", 12); | ||
for (unsigned int i_bc = 0; | ||
i_bc < this->simulation_parameters.boundary_conditions.size; | ||
++i_bc) | ||
{ | ||
this->tracer_flow_rate_table.add_value("flow-rate-" + | ||
Utilities::int_to_string(i_bc, | ||
2), | ||
tracer_flow_rate_vector[i_bc]); | ||
this->tracer_flow_rate_table.set_scientific( | ||
"flow-rate-" + Utilities::int_to_string(i_bc, 2), true); | ||
} | ||
|
||
// Console output | ||
if (simulation_parameters.post_processing.verbosity == | ||
Parameters::Verbosity::verbose) | ||
{ | ||
this->pcout << "Tracer flow rate at the boundaries: " << std::endl; | ||
for (unsigned int i_bc = 0; | ||
i_bc < this->simulation_parameters.boundary_conditions.size; | ||
++i_bc) | ||
this->pcout << "\t boundary " << i_bc << ": " | ||
<< tracer_flow_rate_vector[i_bc] << std::endl; | ||
} | ||
|
||
if ((simulation_control->get_step_number() % | ||
this->simulation_parameters.post_processing.output_frequency == | ||
0) && | ||
Utilities::MPI::this_mpi_process(mpi_communicator) == 0) | ||
{ | ||
std::string filename = | ||
simulation_parameters.simulation_control.output_folder + | ||
simulation_parameters.post_processing.tracer_flow_rate_output_name + | ||
".dat"; | ||
std::ofstream output(filename.c_str()); | ||
this->tracer_flow_rate_table.write_text(output); | ||
} |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Personally I would fill the table write the output outside of the function and not inside. The function calculates the flow rate, it does not really take care of anything else. This should be two functions. One two calculate and one to store the information.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I split into 2 functions as suggested. It makes more sense that way
1fa519b
to
9d076f6
Compare
@oguevremont can you rebase and then I'll be able to merge |
Co-authored-by: Amishga Alphonius <107414376+AmishgaAlphonius@users.noreply.github.com>
Co-authored-by: hepap <47506601+hepap@users.noreply.github.com>
Co-authored-by: Amishga Alphonius <107414376+AmishgaAlphonius@users.noreply.github.com>
9d076f6
to
28693db
Compare
Description An application of the tracer physics is to obtain the residence time distribution (RTD) of a flow through a domain. This PR adds the post-processing option that outputs the tracer flow rate at every boundary, providing the data for RTD analyses. In addition, a bug where the conditions to allow boundary conditions to update their values were too strict was rectified. Each physics is now responsible for deciding when to update their boundary conditions, and this condition was already based on whether boundary conditions were time-dependent. Co-authored-by: Amishga Alphonius <107414376+AmishgaAlphonius@users.noreply.github.com> Co-authored-by: hepap <47506601+hepap@users.noreply.github.com> Former-commit-id: fd4481e
Description An application of the tracer physics is to obtain the residence time distribution (RTD) of a flow through a domain. This PR adds the post-processing option that outputs the tracer flow rate at every boundary, providing the data for RTD analyses. In addition, a bug where the conditions to allow boundary conditions to update their values were too strict was rectified. Each physics is now responsible for deciding when to update their boundary conditions, and this condition was already based on whether boundary conditions were time-dependent. Co-authored-by: Amishga Alphonius <107414376+AmishgaAlphonius@users.noreply.github.com> Co-authored-by: hepap <47506601+hepap@users.noreply.github.com> Former-commit-id: fd4481e
Description
An application of the tracer physics is to obtain the residence time distribution (RTD) of a flow through a domain. This PR adds the post-processing option that outputs the tracer flow rate at every boundary, providing the data for RTD analyses.
In addition, a bug where the conditions to allow
boundary conditions
to update their values were too strict was rectified. Each physics is now responsible for deciding when to update their boundary conditions, and this condition was already based on whether boundary conditions were time-dependent.Testing
Documentation
Miscellaneous (will be removed when merged)
Checklist (will be removed when merged)
See this page for more information about the pull request process.
Code related list:
Pull request related list: