-
Notifications
You must be signed in to change notification settings - Fork 46
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
Saving specific variables in OpenQBMM in OpenFOAM to enable its use later within OpenFOAM #166
Comments
Hi @kurianjv and welcome! If I understand your question correctly, you want the entire history of the EQMOM reconstruction (primary and secondary weights, abscissae, sigma). If that is the case, see my comments below.
This approach only works for post-processing or for restart purposes. If you need to save all the time steps, it will require to change controlDict so that all the fields that have AUTO_WRITE active to be written, which is going to be inefficient and take a lot of space.
This is not very different from the previous approach. It will require to store all the moments at each time step, and then for each step, invoke the moment inversion to reconstruct the NDF with EQMOM again.
In OpenFOAM you can force a field to be written at any point in the code with the "write()" method called on the field. For example:
This will create a directory with the time step value as name, and inside it will only store the fields you forced to write. The downside of this approach is that is it not clean: you still get some directories with all the fields that are saved following controlDict. OpenFOAM will not complain (except if you use the incomplete directory as restart point), but Paraview may be less happy. An alternative for fields registered in the OpenFOAM database can be to use the
I am not sure about the purpose of this, but if you need to compute source terms based on the EQMOM reconstruction, it may easier on the long run to couple the codes. If you provide some more details on what is the final goal, I can try to guide. I hope this helps! |
Hey Alberto Thank you for the quick reply! For example, I was trying to print moment once the populationBalance->solve() is excecuted by I am thinking out loud:) |
Let me try to explain the structure of the code a little bit. First, quadratureNode is a templated class, so you cannot directly use it as you tried to because you are not specifying the type of field for weights and abscissae. For example, This said, if I understand correctly, you just want to access weights and abscissae that are stored already in memory. The There is no direct access function to this class from
To make pbeTransportFoam build with these additions, we need to add
to and modify the
The access functions to weights and abscissae are declared here https://github.com/OpenQBMM/OpenQBMM/blob/master/src/quadratureMethods/quadratureNode/quadratureNode/quadratureNodeI.H Specifically:
Notes
|
Hey Alberto Perfect! That worked :) |
You are welcome. I close based on the comment above. Please, reopen at need. |
Hey Alberto I was trying to access the source terms used for moment transport equation.
I was trying to save these source terms similarly to volScalarFields in pbetransportFoam (like I was trying to save weights and abscissae which was discussed previously). Could you give some thoughts on how I can access these source terms once they are computed to update the each moment? |
Hi @kurianjv the source terms are computed on a cell basis, not as fields, because they are solved with an adaptive ODE solver in each cell. The evaluation of the source terms of univariatePopulationBalance is done in Line 151 in 3f0b4c5
However, there is no mechanism in place to access the field of source terms. Why do you need to access it: to change it or to use the value for processing? |
Hey Alberto I wanted to get the field for processing, if it is not accessible I will add code to compute it accordingly then. Thank you :) where n is defined as n = node.n(celli, primaryweight, bAbscissa) and d = node.d(celli, bAbscissa). Could you help me in understanding what are n and d ? Is this understanding correct ? Thank you so much in advance for taking the time! |
Hi Kurian,
OK. Just one note: if the source you need is the one for aggregation, recomputing it is very expensive, particularly with EQMOM, compared to the rest of the calculation, because of the number of operations involved. What you could do define extra fields in the OpenQBMM code, with a name, and then access to them using the database, as discussed above. Note that the fields should be properties of some object that persists (like the model itself). You could save one cell value at a time in these fields, and then use the fields later.
n is the number density at that cell:
d is the diameter at a cell given the abscissa:
Correct. The pos function is used to define a step function (pos returns a field which is 1 when the argument is positive, zero otherwise). |
Hey again Alberto Thank you for the information about files in which both n and d are defined ! But in Do let me know I have misunderstood something ?
Thanks again for taking the time! |
Hi Kurian, you are correct, b is calculated
in different ways depending on how the NDF is defined and what the weights represents. If the distribution is a number density, then the weight is also a number density (number/m^3). However, in some applications this leads to very large numbers, so it is convenient to redefine the equations so that the weight is a volume fraction. This is typically the case when we deal with multiphase flows with large number of bubbles or particles. This is why the n method takes care of returning the actual number density: because the weight may be a number density or a volume fraction or a mass density, depending on the The growth kernel expects a weight to be a number density, and the abscissa to be a length, that's why the n method is called, so the same kernels can be reused. Note that if volumeFraction = false, then the weight is directly returned:
About the references, the main ones are listed here: https://www.openqbmm.org/references-and-publications/ In summary:
For EQMOM (note these two papers were not produced with OpenQBMM: I clarify because someone in a reproducibility study thought they were for some reason. It is OpenQBMM to be based on them):
The two main publications on OpenQBMM are (preprints are on ResearchGate):
A new one is currently under review, on joint size-velocity NDFs. The preprint is here: https://www.researchgate.net/publication/348538579_A_quadrature-based_moment_method_for_the_evolution_of_the_joint_size-velocity_number_density_function_of_a_particle_population |
Hey Alberto Thank you for the detailed explanations, i understand much better now! |
You are welcome. Good luck with your project! |
Hey
I was hoping to get some help with respect to saving a variable which is computed in OpenQBMM, like the secondary weights and abscissae at every time step. When I looked at the OpenQBMM 5.0.1 and the issues previously raised on this platform, I found two ways to print variables:
I cannot figure out how either of these methods can be customized to allow for direct saving variables in OpenFOAM. My thought was to save each of these variables (in a generic way) like under a volScalarField and then use it later in the OpenFOAM solver (rather than within OpenQBMM). This would allow for me to compute NDF at every point in the domain at every time step or calculate other source terms.
Any help would be greatly appreciated!
And thanks in advance
The text was updated successfully, but these errors were encountered: