OpenTidalFarm is a layout optimisation software for tidal turbine farms.
The positioning of the turbines in a tidal farm is a crucial decision. Simulations show that the optimal positioning can increase the power generation of the farm by up to 50% and can therefore determine the viability of a project.
However, finding the optimal layout is a difficult process due to the complex flow interactions. OpenTidalFarm solves this problem by applying an efficient optimisation algorithm onto a accurate flow prediction model.
Following presentation gives a quick introduction to OpenTidalFarm:
The following video demonstrates OpenTidalFarm on different setups:
The following videos show the how the turbine positions change during the optimisation. The final frame of the videos show the optimisation locations of the turbines. In these cases, the total power production of the turbines was improved by over 25%.
- Optimisation of 128 turbines in the Pentland Firth, Scotland
- Optimisation of 128 turbines in the Pentland Firth, Scotland (zoomed into site area)
- Optimisation of 256 turbines in the Pentland Firth, Scotland
- Optimisation of 256 turbines in the Pentland Firth, Scotland (zoomed into site area)
- High resolution shallow water model for accurate flow prediction.
- Arbitrary shoreline data and bathymetry support.
- Optimise the turbine position and size to maximise the total farm power output.
- Site constraints / minimum distance between turbines.
- Optimisation of up to hundreds of turbines.
- Checkpoint support to restart optimisation.
Note: This installation procedure assumes that you are running the Ubuntu operating system.
The installation consists of following steps
- Download and install the dependencies:
- FEniCS project >=1.2 (Follow the Ubuntu PPA installation order to get a recent installation)
- dolfin-adjoint.
- SciPy Version >=0.11 - e.g. with
pip install scipy
.
- Download OpenTidalFarm and extract it.
- Open a terminal and change into the extracted directory and run
sudo python setup.py install
to install it. A simple test to check if the installation was correct is to open a Python shell and type:
from opentidalfarm import *
If you get an error, make sure that you have set the PYTHONPATH
correctly. In Linux, this can be done with:
export PYTHONPATH=/XYZ:$PYTHONPATH
where XYZ
should be replaced with the path to your OpenTidalFarm installation.
Now you are ready to run one of the many examples in the examples/
folder.
Following example code shows how to optimise the position of 32 turbines in a mesh of the Orkney islands.
from opentidalfarm import *
config = SteadyConfiguration("mesh/earth_orkney_converted.xml", inflow_direction=[0.9865837220518425, -0.16325611591095968])
config.params['diffusion_coef'] = 90.0
config.params['turbine_x'] = 40.
config.params['turbine_y'] = 40.
config.params['controls'] = ['turbine_pos']
# Some domain information extracted from the geo file.
# This information is used to deploy the turbines autmatically.
site_x = 1000.
site_y = 500.
site_x_start = 1.03068e+07
site_y_start = 6.52246e+06 - site_y
config.set_site_dimensions(site_x_start, site_x_start + site_x, site_y_start, site_y_start + site_y)
# Place 32 turbines in a regular grid, each with a maximum friction coefficient of 10.5
deploy_turbines(config, nx=8, ny=4, friction=10.5)
# Define some constraints for the optimisation positions.
# Constraint to keep the turbines within the site area
lb, ub = position_constraints(config)
# Constraint to keep a minium distance of 1.5 turbine diameter between each turbine
ineq = get_minimum_distance_constraint_func(config)
# Solve the optimisation problem
rf = ReducedFunctional(config, plot=True)
maximize(rf, bounds=[lb, ub], constraints=ineq, method="SLSQP")
This example can be found in the examples/tutorial
directory and can be executed by running make mesh && make
.
The output files are:
- turbine.pvd: The turbine positions at each optimisation step
- p2p1_u.pvd: The velocity function for the most recent turbine position calculation.
- p2p1_p.pvd: The free-surface displacement function for the most recent turbine position calculation.
If you only want to compute the power production for the given layout (without optimising), replace the ast code line above with:
# Switch off the computation of the automatic scaling factor (requires one adjoint solve), as it is needed only for the optimisation
config.params['automatic_scaling'] = False
# Retrieve the initial control values (here: turbine positions)
m0 = rf.initial_control()
# Store the turbine positions as a pvd file
rf.update_turbine_cache(m0)
File("turbines.pvd") << config.turbine_cache.cache["turbine_field"]
# Compute the extracted power from the flow
print "Power extraction: %f MW" % rf.j(m0)
OpenTidalFarm is based on configurations for defining different setups. The most important configurations are
SteadyConfiguration
: Use this configuration for steady simulations.UnstadyConfiguration
: Use this configuration for unsteady simulations.
Once you have a configuration object, try running
config.info()
to see detailed information about the settings of the current configuration.
Each configuration has a large of additional parameters that can be changed.
For example to use Picard iterations instead of a Newton solver to solve the nonlinear shallow water equations one would set:
config.params["newton_solver"] = False
Some of the more important parameters are:
- "controls": Defines the control parameters that the optimisation algorithm may use. Possible choicees are the optimisation of the turbine positions and/or the friction of each individual turbine. Valid values: a list containing one or more of
['turbine_pos', 'turbine_friction']
. - "save_checkpoints": Automatically save checkpoints to disk from which the optimisation can be restarted. Valid values:
True
orFalse
. - "turbine_x": The x-extension of each turbine.
- "turbine_y": The y-extension of each turbine.
- "output_turbine_power": Output the energy production for each individual turbine. Valid values:
True
orFalse
.
Again, use config.info()
to list the current configuration setup.
Often, one ones to ensure that the turbines must be deployed in a particular area. This can be achieved in two ways:
Box constraints is the easiest solution of the deployment area is a rectangle.
The code below shows an example usage:
config.set_site_dimensions(site_x_start, site_x_end, site_y_start, site_y_end)
bounds = position_constraints(config)
maximize(rf, bounds=bounds, method = "SLSQP")
For more complex site domains, one can define a polygon in which the turbines must be deployed.
The code below shows an example usage for a domain polygon with 4 vertices:
vertices = [[site_x_start, site_y_start],
[site_x_end, site_y_start],
[site_x_end, site_y_end],
[site_x_start, site_y_end]]
ineq = generate_site_constraints(config, vertices)
maximize(rf, constraints=ineq, method = "SLSQP")
Note that SLSQP in an infeasible algorithm, which means that the initial turbine layout does not have to fulfill these constraints.
OpenTidalFarm expects the 3 the mesh to have three identifiers of the boundary mesh:
- ID 1: inflow boundary
- ID 2: outflow boundary
- ID 3: shoreline boundary
The available optimisation options are explained in detail here.
OpenTidalFarm can automatically store checkpoints to disk from which the optimisation procedure can be restarted. The checkpoints generation is activated with:
config.params["save_checkpoints"] = True
where config
is the Configuration
object.
The checkpoint data is stored in the two files "checkpoint_fwd.dat" and "checkpoint_adj.dat".
In order to restart from a checkpoint you need to load in the checkpoint with:
rf.load_checkpoint()
where rf
is the ReducedFunctionalObject
.
You will see that the optimisation starts from the beginning, however the optimisation iterations until the checkpoint will happen instantly since the solutions are cached in the checkpoint.
By default, OpenTidalFarm only uses the -O3
compiler optimisation flag as a safe choice.
However, for large optimisations runs, one might want to use more aggressive compiler optimisation. Experiments have shown that following options can yield significant speed up:
parameters['form_compiler']['cpp_optimize_flags'] = '-O3 -ffast-math -march=native'
or the less aggressive version:
parameters['form_compiler']['cpp_optimize_flags'] = '-O3 -fno-math-errno -march=native'
Add this line just before you call the maximize function.
However, note that in some circumstances such aggressive optimisation might be problematic for the optimisation algorithms. If the optimisation algorithm returns errors saying that the gradient
Sometimes, the SLSQP optimisation algorithm stops early giving the error:
Positive directional derivative for linesearch (Exit mode 8)
In such case you can try following things:
- If you have the
-ffast-math -march=native
compiler flags active (see above), try switching them off. - Use finer mesh in the turbine site area. The numerical errors of representing the turbines might be dominating the problem.
- Use a looser optimisation tolerance, by passing the
tol
parameter to maximize function. - If you are only optimising for the turbine friction and you do not use
automatic_scaling
parameter (which you should in this particular case), then you might have to scale your problem manually. This is done by passing a scale argument to maximize, e.g.maximize(rf, scale=1e-3)
. Choose the scale value such that the first gradient is in the order of 1.
OpenTidalFarm updates the configuration values while running the optimisation. Therefore, in order to get the position and friction values after the optimisation, simply add this line of code at the bottom of your program:
print "Turbine positions: ", config.params['turbine_pos']
print "Turbine friction: ", config.params['turbine_friction']
The location, power generated and turbine friction of the optimised layout can be stored in .csv files named 'turbine_info.csv' in the iteration directories, simply add this line to your program:
config.params["print_individual_turbine_power"] = True
Additionally, the very results of the initial guess is stored in 'iter_0/initial_turbine_info.csv'.
Please cite the following paper if you are using OpenTidalFarm:
S.W. Funke, P.E. Farrell, M.D. Piggott, Tidal turbine array optimisation using the adjoint approach, Renewable Energy, accepted (2013) arXiv:1304.1768
For the automated optimisation framework used by OpenTidalFarm, please cite:
Simon W. Funke and Patrick E. Farrell. A framework for automated PDE-constrained optimisation, TOMS, submitted. arXiv:1302.3894
For the automatic adjoint generation used by OpenTidalFarm, please cite:
Patrick E. Farrell, David A. Ham, Simon W. Funke and Marie E. Rognes (2013). Automated derivation of the adjoint of high-level transient finite element programs. SIAM Journal on Scientific Computing, Vol:35, ISSN:1064-8275, Pages:C369-C393
For questions and support please contact Simon Funke [email protected].
OpenTidalFarm is an open source project that can be freely used under the GNU GPL version 3 licence.