Calculation of thermodynamic properties and phase-equilibria with Equation of State. Keeping high performance with the power of Fortran and the possiblity of using either automatic differentiation and anallytical derivatives.
Find us on…
There are multiple open source equations of state libraries, like:
julia
rust
with Python
bindingsC++
with Python
bindingspython
Fortran
with Python
bindingsC++
with Python
bindingsHere we are presenting yet another (still in development) one, that tackles the same problem just, in another way. Mostly exploiting the readability and extensibility of Modern Fortran for scientists to have an easy way to implement new thermodynamic models without dealing with lower-level languages but still getting decent performance. And also this framework provides the possibility of using analytically obtained derivatives, so both options are easily available.
This is an experimental work in progress, and we recommend the before mentioned libraries if you are intending to use some of this in real work. Big part of the code comes from a refactoring process of older codes so not all parts are easily readable, yet.
We focus mainly on that the addition of a new thermodynamic model as easily as possible. Also providing our models too!
The latest API documentation for the main
branch can be found:
- Fortran API documentation
- Python API documentation
The Fortran API documentation can also be generated by processing the source files with FORD. On the other hand, the Python API documentation can be generated by processing the source files with Sphinx.
This library is currently maintained by the research group of Prof. MartÃn Cismondi-Duarte at IPQA (UNC-CONICET)
yaeos
A lot of users get the bad picture of Fortran being old and archaic since most
of the codes they’ve seen are written in ancient F77
.
use yaeos, only: PengRobinson76, ArModel
integer, parameter :: n=2 ! Number of components
real(8) :: V, T, P, dPdN(n) ! variables to calculate
class(ArModel), allocatable :: model ! Model
real(pr) :: z(n), tc(n), pc(n), w(n), kij(n, n), lij(n, n)
z = [0.3, 0.7]
tc = [190., 310.]
pc = [14., 30.]
w = [0.001, 0.03]
kij = reshape([0., 0.1, 0.1, 0.], [n,n])
lij = kij / 2
model = PengRobinson76(tc, pc, w, kij, lij)
V = 1
T = 150
call model%pressure(z, V, T, P)
print *, P
! Obtain derivatives adding them as optional arguments!
call model%pressure(model, z, V, T, P, dPdN=dPdN)
print *, dPdN
Examples of code with simple applications showing the capabilities of yaeos
can be found at example/tutorials. Each example can be run
with:
fpm run --example <example name here>
Not providing any examples will show all the possible examples that can be run.
yaeos
needs to have both lapack
and nlopt
libraries on your system.
sudo apt install libnlopt-dev libblas-dev liblapack-dev
yaeos
yaeos
is intended to use as a fpm
package. fpm
is the Fortran Package Manager, which automates the compilation and running
process of Fortran libraries/programs.
You can either:
yaeos
as a dependency with:
bash fpm new my_project
In the
fpm.toml
file add:
toml [dependencies] yaeos = {git="https://github.com/ipqa-research/yaeos"}
app
directory
bash git clone https://github.com/ipqa-research/yaeos cd yaeos fpm run
If your intention is either to develop for yaeos
or to explore in more detail
the library with debugging. We provide some predefined defuaults to work with
vscode
. You can add them to the cloned repository by running:
git clone https://github.com/ipqa-research/vscode-fortran .vscode
From the project main directory
In this repository we provide a series of examples of the different things that
can be calculated with yaeos
. The source codes for the examples can be seen
at the example/tutorials directory.
All the examples can be run with
fpm run --example <example_name_here>
We are using the hyperdual
module developed by
Philipp Rehner
and Gernot Bauer
The automatic differentiation API isn’t fully optimized yet so performance is much slower than it should be.
A complete implementation of the PR76 Equation of State can me found in
example/adiff/adiff_pr76.f90
. Or in the documentation pages.
It is also possible to differentiate with tapenade
. Examples can be seen
in the documentation pages or in The tools directory