From 8181b97cd1e3cc1e783ecc262ab19ab929f148be Mon Sep 17 00:00:00 2001 From: richardreeve Date: Mon, 15 Apr 2024 00:30:48 +0100 Subject: [PATCH 1/3] Move MPI structs into extension --- src/EcoSISTEM.jl | 41 +++--------------- src/MPIEcosystem.jl | 103 +++++++++++++++++++++++++++++--------------- src/MPIGenerate.jl | 86 ++++++++++++++++++++---------------- src/MPILandscape.jl | 90 ++++++++++++++++++++++---------------- 4 files changed, 176 insertions(+), 144 deletions(-) diff --git a/src/EcoSISTEM.jl b/src/EcoSISTEM.jl index 05557aeb..9a74ea26 100644 --- a/src/EcoSISTEM.jl +++ b/src/EcoSISTEM.jl @@ -94,42 +94,11 @@ export meta_simpson, meta_shannon, meta_speciesrichness, mean_abun, using Random -""" - MPIGridLandscape{RA <: Base.ReshapedArray, NT <: NamedTuple} - -MPIEcosystem abundances housed in the landscape, shared across multiple nodes. -""" -mutable struct MPIGridLandscape{RA <: Base.ReshapedArray, NT <: NamedTuple} - rows_matrix::Matrix{Int64} - cols_vector::Vector{Int64} - reshaped_cols::Vector{RA} - rows_tuple::NT - cols_tuple::NT - rngs::Vector{MersenneTwister} -end - -""" - MPIEcosystem{MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, - SL <: SpeciesList, TR <: AbstractTraitRelationship} <: AbstractEcosystem{Part, SL, TR} - -MPIEcosystem houses information on species and their interaction with their environment. It houses all information of a normal `Ecosystem` (see documentation for more details), with additional fields to describe which species are calculated on which machine. This includes: `sppcounts` - a vector of number of species per node, `firstsp` - the identity of the first species held by that particular node. -""" -mutable struct MPIEcosystem{MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, - SL <: SpeciesList, - TR <: AbstractTraitRelationship} <: - AbstractEcosystem{Part, SL, TR} - abundances::MPIGL - spplist::SL - abenv::Part - ordinariness::Union{Matrix{Float64}, Missing} - relationship::TR - lookup::Vector{Lookup} - sppcounts::Vector{Int32} - firstsp::Int64 - sccounts::Vector{Int32} - firstsc::Int64 - cache::Cache -end +abstract type MPIGridLandscape end +abstract type MPIEcosystem{MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, + SL <: SpeciesList, + TR <: AbstractTraitRelationship} <: + AbstractEcosystem{Part, SL, TR} end export MPIGridLandscape export MPIEcosystem diff --git a/src/MPIEcosystem.jl b/src/MPIEcosystem.jl index 94ddab6c..3abaa0c8 100644 --- a/src/MPIEcosystem.jl +++ b/src/MPIEcosystem.jl @@ -1,4 +1,4 @@ -using EcoSISTEM +import EcoSISTEM using MPI using Diversity using HCubature @@ -6,40 +6,71 @@ using Unitful using EcoSISTEM.Units using Missings -import EcoSISTEM: MPIEcosystem, gather_abundance, gather_diversity -using EcoSISTEM: AbstractAbiotic, AbstractTraitRelationship, Lookup, Cache - -function MPIEcosystem(abundances::MPIGL, - spplist::SL, abenv::Part, - ordinariness::Union{Matrix{Float64}, Missing}, - relationship::TR, lookup::Vector{Lookup}, - sppcounts::Vector, firstsp::Int64, sccounts::Vector, - firstsc::Int64, - cache::Cache) where - {MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, - SL <: SpeciesList, TR <: AbstractTraitRelationship} - tematch(spplist, abenv) || error("Traits do not match habitats") - trmatch(spplist, relationship) || - error("Traits do not match trait functions") - #_mcmatch(abundances.matrix, spplist, abenv) || - # error("Dimension mismatch") - return MPIEcosystem{MPIGL, Part, SL, TR}(abundances, spplist, abenv, - ordinariness, - relationship, lookup, sppcounts, - firstsp, - sccounts, firstsc, cache) +""" + MPIEcosystem{MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, + SL <: SpeciesList, TR <: AbstractTraitRelationship} <: + AbstractEcosystem{Part, SL, TR} + +MPIEcosystem houses information on species and their interaction with their +environment. It houses all information of a normal `Ecosystem` (see +documentation for more details), with additional fields to describe which +species are calculated on which machine. This includes: `sppcounts` - a +vector of number of species per node, `firstsp` - the identity of the +first species held by that particular node. +""" +mutable struct MPIEcosystem{MPIGL <: EcoSISTEM.MPIGridLandscape, + Part <: EcoSISTEM.AbstractAbiotic, + SL <: EcoSISTEM.SpeciesList, + TR <: EcoSISTEM.AbstractTraitRelationship} <: + EcoSISTEM.MPIEcosystem{MPIGL, Part, SL, TR} + abundances::MPIGL + spplist::SL + abenv::Part + ordinariness::Union{Matrix{Float64}, Missing} + relationship::TR + lookup::Vector{EcoSISTEM.Lookup} + sppcounts::Vector{Int32} + firstsp::Int64 + sccounts::Vector{Int32} + firstsc::Int64 + cache::EcoSISTEM.Cache + + function MPIEcosystem(abundances::MPIGL, + spplist::SL, abenv::Part, + ordinariness::Union{Matrix{Float64}, + Missing}, + relationship::TR, + lookup::Vector{EcoSISTEM.Lookup}, + sppcounts::Vector, + firstsp::Int64, + sccounts::Vector, + firstsc::Int64, + cache::EcoSISTEM.Cache) where {MPIGL, Part, SL, TR} + EcoSISTEM.tematch(spplist, abenv) || error("Traits do not match habitats") + EcoSISTEM.trmatch(spplist, relationship) || + error("Traits do not match trait functions") + #_mcmatch(abundances.matrix, spplist, abenv) || + # error("Dimension mismatch") + return new{MPIGL, Part, SL, TR}(abundances, spplist, abenv, + ordinariness, relationship, lookup, + sppcounts, firstsp, sccounts, firstsc, + cache) + end end +EcoSISTEM.MPIEcosystem(args...) = MPIEcosystem(args...) + using EcoSISTEM: getkernels, genlookups, numrequirements """ MPIEcosystem(spplist::SpeciesList, abenv::GridAbioticEnv, rel::AbstractTraitRelationship) -Function to create an `MPIEcosystem` given a species list, an abiotic environment and trait relationship. +Function to create an `MPIEcosystem` given a species list, an abiotic +environment and trait relationship. """ -function MPIEcosystem(popfun::F, spplist::SpeciesList{T, Req}, - abenv::GridAbioticEnv, - rel::AbstractTraitRelationship) where +function MPIEcosystem(popfun::F, spplist::EcoSISTEM.SpeciesList{T, Req}, + abenv::EcoSISTEM.GridAbioticEnv, + rel) where {F <: Function, T, Req} comm = MPI.COMM_WORLD rank = MPI.Comm_rank(comm) @@ -60,7 +91,7 @@ function MPIEcosystem(popfun::F, spplist::SpeciesList{T, Req}, firstsc = scindices[rank + 1] + 1 # Create matrix landscape of zero abundances - ml = emptyMPIgridlandscape(sppcounts, sccounts) + ml = EcoSISTEM.emptyMPIgridlandscape(sppcounts, sccounts) # Populate this matrix with species abundances popfun(ml, spplist, abenv, rel) @@ -71,12 +102,14 @@ function MPIEcosystem(popfun::F, spplist::SpeciesList{T, Req}, nm = zeros(Int64, (sppcounts[rank + 1], numsc)) totalE = zeros(Float64, (numsc, numrequirements(Req))) return MPIEcosystem(ml, spplist, abenv, missing, rel, lookup_tab, sppcounts, - firstsp, sccounts, firstsc, Cache(nm, totalE, false)) + firstsp, sccounts, firstsc, + EcoSISTEM.Cache(nm, totalE, false)) end -function MPIEcosystem(spplist::SpeciesList, abenv::GridAbioticEnv, - rel::AbstractTraitRelationship) - return MPIEcosystem(populate!, spplist, abenv, rel) +function MPIEcosystem(spplist::EcoSISTEM.SpeciesList, + abenv::EcoSISTEM.GridAbioticEnv, + rel) + return MPIEcosystem(EcoSISTEM.populate!, spplist, abenv, rel) end """ @@ -84,7 +117,7 @@ end Gather full abundances matrix on root node. """ -function gather_abundance(eco::MPIEcosystem) +function EcoSISTEM.gather_abundance(eco::MPIEcosystem) comm = MPI.COMM_WORLD rank = MPI.Comm_rank(comm) true_abuns = zeros(Int64, counttypes(eco), countsubcommunities(eco)) @@ -143,8 +176,8 @@ function _calcordinariness(eco::MPIEcosystem) relab end -function gather_diversity(eco::MPIEcosystem, divmeasure::F, - q) where {F <: Function} +function EcoSISTEM.gather_diversity(eco::MPIEcosystem, divmeasure::F, + q) where {F <: Function} comm = MPI.COMM_WORLD rank = MPI.Comm_rank(comm) totalsize = MPI.Comm_size(comm) diff --git a/src/MPIGenerate.jl b/src/MPIGenerate.jl index 93a64640..4b303b6a 100644 --- a/src/MPIGenerate.jl +++ b/src/MPIGenerate.jl @@ -1,19 +1,18 @@ -using EcoSISTEM +import EcoSISTEM using MPI using LinearAlgebra using Distributions -import EcoSISTEM: getlookup, update!, update_energy_usage!, move!, populate! using EcoSISTEM: AbstractAbiotic, Abstract1Requirement, Abstract2Requirements using EcoSISTEM: AbstractHabitat, AbstractBudget, AbstractTraitRelationship using EcoSISTEM: energy_adjustment, invalidatecaches!, habitatupdate!, - budgetupdate! + budgetupdate!, BirthOnlyMovement, BudgetCollection2 """ update!(eco::MPIEcosystem, timestep::Unitful.Time) where N Function to update an MPIEcosystem abundances and environment for one timestep. """ -function update!(eco::MPIEcosystem, timestep::Unitful.Time) +function EcoSISTEM.update!(eco::MPIEcosystem, timestep::Unitful.Time) comm = MPI.COMM_WORLD rank = MPI.Comm_rank(comm) @@ -21,7 +20,7 @@ function update!(eco::MPIEcosystem, timestep::Unitful.Time) numsc = countsubcommunities(eco) params = eco.spplist.params # Set the overall energy budget of that square - update_energy_usage!(eco) + EcoSISTEM.update_energy_usage!(eco) MPI.Allgatherv!(MPI.VBuffer(eco.cache.totalE, eco.sccounts), comm) eco.cache.valid = true @@ -37,7 +36,7 @@ function update!(eco::MPIEcosystem, timestep::Unitful.Time) sc, truesp) # Convert 1D dimension to 2D coordinates - (x, y) = convert_coords(eco, sc) + (x, y) = EcoSISTEM.convert_coords(eco, sc) # Check if grid cell currently active if eco.abenv.active[x, y] && (eco.cache.totalE[sc, 1] > 0) @@ -63,7 +62,7 @@ function update!(eco::MPIEcosystem, timestep::Unitful.Time) eco.abundances.rows_matrix[mpisp, sc] += (births - deaths) # Calculate moves and write to cache - move!(eco, eco.spplist.movement, sc, truesp, + EcoSISTEM.move!(eco, eco.spplist.movement, sc, truesp, eco.cache.netmigration, births) end end @@ -71,7 +70,7 @@ function update!(eco::MPIEcosystem, timestep::Unitful.Time) # Update abundances with all movements eco.abundances.rows_matrix .+= eco.cache.netmigration - synchronise_from_rows!(eco.abundances) + EcoSISTEM.synchronise_from_rows!(eco.abundances) # Invalidate all caches for next update invalidatecaches!(eco) @@ -81,13 +80,17 @@ function update!(eco::MPIEcosystem, timestep::Unitful.Time) return budgetupdate!(eco, timestep) end -function getlookup(eco::MPIEcosystem, sp::Int64) +function EcoSISTEM.getlookup(eco::MPIEcosystem, sp::Int64) return eco.lookup[sp - eco.firstsp + 1] end -function update_energy_usage!(eco::MPIEcosystem{MPIGL, A, - SpeciesList{Tr, Req, B, C, D}, - E}) where +function EcoSISTEM.update_energy_usage!(eco::MPIEcosystem{MPIGL, A, + EcoSISTEM.SpeciesList{Tr, + Req, + B, + C, + D}, + E}) where {MPIGL <: MPIGridLandscape, A, B, C, D, E, Tr, Req <: Abstract1Requirement} !eco.cache.valid || return true @@ -115,15 +118,21 @@ function update_energy_usage!(eco::MPIEcosystem{MPIGL, A, return eco.cache.valid = true end -function update_energy_usage!(eco::MPIEcosystem{MPIGL, A, - SpeciesList{Tr, Req, B, C, D}, - E}) where { - MPIGL <: - MPIGridLandscape, A, - B, C, D, E, Tr, - Req <: - Abstract2Requirements - } +function EcoSISTEM.update_energy_usage!(eco::MPIEcosystem{MPIGL, A, + EcoSISTEM.SpeciesList{Tr, + Req, + B, + C, + D}, + E}) where { + MPIGL <: + MPIGridLandscape, + A, + B, C, D, E, + Tr, + Req <: + Abstract2Requirements + } !eco.cache.valid || return true rank = MPI.Comm_rank(MPI.COMM_WORLD) @@ -155,11 +164,12 @@ function update_energy_usage!(eco::MPIEcosystem{MPIGL, A, end using EcoSISTEM: getdimension, getboundary, calc_lookup_moves! -function move!(eco::MPIEcosystem, ::BirthOnlyMovement, sc::Int64, truesp::Int64, - grd::Matrix{Int64}, births::Int64) +function EcoSISTEM.move!(eco::MPIEcosystem, ::BirthOnlyMovement, sc::Int64, + truesp::Int64, + grd::Matrix{Int64}, births::Int64) width, height = getdimension(eco) - (x, y) = convert_coords(eco, sc, width) - lookup = getlookup(eco, truesp) + (x, y) = EcoSISTEM.convert_coords(eco, sc, width) + lookup = EcoSISTEM.getlookup(eco, truesp) calc_lookup_moves!(getboundary(eco.spplist.movement), x, y, truesp, eco, births) # Lose moves from current grid square @@ -170,16 +180,17 @@ function move!(eco::MPIEcosystem, ::BirthOnlyMovement, sc::Int64, truesp::Int64, for i in eachindex(lookup.x) newx = mod(lookup.x[i] + x - 1, width) + 1 newy = mod(lookup.y[i] + y - 1, height) + 1 - loc = convert_coords(eco, (newx, newy), width) + loc = EcoSISTEM.convert_coords(eco, (newx, newy), width) grd[mpisp, loc] += mov[i] end return eco end using EcoSISTEM: _getdimension, _getbudget -function populate!(ml::MPIGridLandscape, spplist::SpeciesList, abenv::AB, - rel::R) where {AB <: AbstractAbiotic, - R <: AbstractTraitRelationship} +function EcoSISTEM.populate!(ml::MPIGridLandscape, + spplist::EcoSISTEM.SpeciesList, abenv::AB, + rel::R) where {AB <: AbstractAbiotic, + R <: AbstractTraitRelationship} dim = _getdimension(abenv.habitat) len = dim[1] * dim[2] grid = collect(1:len) @@ -195,14 +206,17 @@ function populate!(ml::MPIGridLandscape, spplist::SpeciesList, abenv::AB, rand!(Multinomial(abundances[mpisp], B), (@view ml.rows_matrix[mpisp, :])) end - return synchronise_from_rows!(ml) + return EcoSISTEM.synchronise_from_rows!(ml) end -function populate!(ml::MPIGridLandscape, spplist::SpeciesList, - abenv::GridAbioticEnv{H, BudgetCollection2{B1, B2}}, - rel::R) where {H <: AbstractHabitat, B1 <: AbstractBudget, - B2 <: AbstractBudget, - R <: AbstractTraitRelationship} +function EcoSISTEM.populate!(ml::MPIGridLandscape, + spplist::EcoSISTEM.SpeciesList, + abenv::EcoSISTEM.GridAbioticEnv{H, + BudgetCollection2{B1, B2}}, + rel::R) where {H <: AbstractHabitat, + B1 <: AbstractBudget, + B2 <: AbstractBudget, + R <: AbstractTraitRelationship} # Calculate size of habitat dim = _getdimension(abenv.habitat) len = dim[1] * dim[2] @@ -222,5 +236,5 @@ function populate!(ml::MPIGridLandscape, spplist::SpeciesList, rand!(Multinomial(abundances[mpisp], B ./ sum(B)), (@view ml.rows_matrix[mpisp, :])) end - return synchronise_from_rows!(ml) + return EcoSISTEM.synchronise_from_rows!(ml) end diff --git a/src/MPILandscape.jl b/src/MPILandscape.jl index f28db3a2..452feb72 100644 --- a/src/MPILandscape.jl +++ b/src/MPILandscape.jl @@ -1,55 +1,71 @@ -using EcoSISTEM +import EcoSISTEM using MPI using Random -import EcoSISTEM: MPIGridLandscape, emptyMPIgridlandscape, - synchronise_from_rows!, synchronise_from_cols! +""" + MPIGridLandscape{RA <: Base.ReshapedArray, NT <: NamedTuple} -function MPIGridLandscape(sppcounts::Vector{Int32}, sccounts::Vector{Int32}, - rows_matrix::Matrix{Int64}, - cols_vector::Vector{Int64}) - rank = MPI.Comm_rank(MPI.COMM_WORLD) +MPIEcosystem abundances housed in the landscape, shared across multiple nodes. +""" +mutable struct MPIGridLandscape{RA <: Base.ReshapedArray, NT <: NamedTuple} <: + EcoSISTEM.MPIGridLandscape + rows_matrix::Matrix{Int64} + cols_vector::Vector{Int64} + reshaped_cols::Vector{RA} + rows_tuple::NT + cols_tuple::NT + rngs::Vector{MersenneTwister} - totalspp = sum(sppcounts) - totalsc = sum(sccounts) + function MPIGridLandscape(sppcounts::Vector{Int32}, + sccounts::Vector{Int32}, + rows_matrix::Matrix{Int64}, + cols_vector::Vector{Int64}) + rank = MPI.Comm_rank(MPI.COMM_WORLD) - lastsp = sum(sppcounts[1:(rank + 1)]) - firstsp = lastsp - sppcounts[rank + 1] + 1 + totalspp = sum(sppcounts) + totalsc = sum(sccounts) - lastsc = sum(sccounts[1:(rank + 1)]) - firstsc = lastsc - sccounts[rank + 1] + 1 + lastsp = sum(sppcounts[1:(rank + 1)]) + firstsp = lastsp - sppcounts[rank + 1] + 1 - sppindices = [0; cumsum(sppcounts) .* sccounts[rank + 1]] - scindices = [0; cumsum(sccounts) .* sppcounts[rank + 1]] + lastsc = sum(sccounts[1:(rank + 1)]) + firstsc = lastsc - sccounts[rank + 1] + 1 - reshaped_cols = map(eachindex(sccounts)) do i - return reshape(view(cols_vector, (sppindices[i] + 1):sppindices[i + 1]), - Int64(sppcounts[i]), Int64(sccounts[rank + 1])) + sppindices = [0; cumsum(sppcounts) .* sccounts[rank + 1]] + scindices = [0; cumsum(sccounts) .* sppcounts[rank + 1]] + + reshaped_cols = map(eachindex(sccounts)) do i + return reshape(view(cols_vector, + (sppindices[i] + 1):sppindices[i + 1]), + Int64(sppcounts[i]), Int64(sccounts[rank + 1])) + end + rows = (total = totalspp, first = firstsp, last = lastsp, + counts = sccounts .* sppcounts[rank + 1]) + cols = (total = totalsc, first = firstsc, last = lastsc, + counts = sppcounts .* sccounts[rank + 1]) + (rows.last - rows.first + 1) * cols.total == length(rows_matrix) || + error("rows_matrix size mismatch: $(rows.last - rows.first + 1) * $(cols.total) !=$(length(rows_matrix))") + (cols.last - cols.first + 1) * rows.total == length(cols_vector) || + error("cols_vector size mismatch: $(cols.last - cols.first + 1) * $(rows.total) !=$(length(cols_vector))") + + return new{typeof(reshaped_cols[1]), typeof(rows)}(rows_matrix, + cols_vector, + reshaped_cols, + rows, cols, + [MersenneTwister(rand(UInt)) + for _ in Base.OneTo(Threads.nthreads())]) end - rows = (total = totalspp, first = firstsp, last = lastsp, - counts = sccounts .* sppcounts[rank + 1]) - cols = (total = totalsc, first = firstsc, last = lastsc, - counts = sppcounts .* sccounts[rank + 1]) - (rows.last - rows.first + 1) * cols.total == length(rows_matrix) || - error("rows_matrix size mismatch: $(rows.last - rows.first + 1) * $(cols.total) !=$(length(rows_matrix))") - (cols.last - cols.first + 1) * rows.total == length(cols_vector) || - error("cols_vector size mismatch: $(cols.last - cols.first + 1) * $(rows.total) !=$(length(cols_vector))") - - return MPIGridLandscape{typeof(reshaped_cols[1]), typeof(rows)}(rows_matrix, - cols_vector, - reshaped_cols, - rows, cols, - [MersenneTwister(rand(UInt)) - for _ in Base.OneTo(Threads.nthreads())]) end +EcoSISTEM.MPIGridLandscape(args...) = MPIGridLandscape(args...) + """ emptyMPIgridlandscape(sppcounts::Vector{Int32}, sccounts::Vector{Int32}) Function to create an empty MPIGridLandscape given information about the MPI setup. """ -function emptyMPIgridlandscape(sppcounts::Vector{Int32}, - sccounts::Vector{Int32}) +function EcoSISTEM.emptyMPIgridlandscape(sppcounts::Vector{Int32}, + sccounts::Vector{Int32}) rank = MPI.Comm_rank(MPI.COMM_WORLD) rows_matrix = zeros(Int64, sppcounts[rank + 1], sum(sccounts)) @@ -58,13 +74,13 @@ function emptyMPIgridlandscape(sppcounts::Vector{Int32}, return MPIGridLandscape(sppcounts, sccounts, rows_matrix, cols_vector) end -function synchronise_from_rows!(ml::MPIGridLandscape) +function EcoSISTEM.synchronise_from_rows!(ml::MPIGridLandscape) return MPI.Alltoallv!(MPI.VBuffer(ml.rows_matrix, ml.rows_tuple.counts), MPI.VBuffer(ml.cols_vector, ml.cols_tuple.counts), MPI.COMM_WORLD) end -function synchronise_from_cols!(ml::MPIGridLandscape) +function EcoSISTEM.synchronise_from_cols!(ml::MPIGridLandscape) return MPI.Alltoallv!(MPI.VBuffer(ml.cols_vector, ml.cols_tuple.counts), MPI.VBuffer(ml.rows_matrix, ml.rows_tuple.counts), MPI.COMM_WORLD) From 4a5169094bf7b6d011c8a375bd2631e7c8c851be Mon Sep 17 00:00:00 2001 From: richardreeve Date: Mon, 15 Apr 2024 00:31:02 +0100 Subject: [PATCH 2/3] Fix organisation name --- CONTRIBUTING.md | 4 ++-- 1 file changed, 2 insertions(+), 2 deletions(-) diff --git a/CONTRIBUTING.md b/CONTRIBUTING.md index 2a68d880..4b0779a6 100644 --- a/CONTRIBUTING.md +++ b/CONTRIBUTING.md @@ -8,9 +8,9 @@ Thanks for contributing to EcoSISTEM.jl! Please read the below on how best to ma 2. Name the branch `username/featurename` and have fun! -3. When your feature is ready to merge back into the `dev` branch create a PR against `boydorr/EcoSISTEM.jl` **and assign `boydorr/EcoSISTEM-admins` as reviewers**. +3. When your feature is ready to merge back into the `dev` branch create a PR against `EcoJulia/EcoSISTEM.jl` **and assign `EcoJulia/EcoSISTEM-admins` as reviewers**. -4. Semver will be handled in PRs from `dev` into `master`. +4. Semver will be handled in PRs from `dev` into `main`. 5. Thanks for your time and effort! From 585b842e156ccbf256455430c534f413d4ab21ef Mon Sep 17 00:00:00 2001 From: richardreeve Date: Mon, 15 Apr 2024 02:46:07 +0100 Subject: [PATCH 3/3] Bump version and news --- NEWS.md | 2 ++ Project.toml | 2 +- src/EcoSISTEM.jl | 6 ++++-- 3 files changed, 7 insertions(+), 3 deletions(-) diff --git a/NEWS.md b/NEWS.md index 4f82315d..128ba44b 100644 --- a/NEWS.md +++ b/NEWS.md @@ -1,5 +1,7 @@ # NEWS +- v0.2.1 + - Move MPI structs into extension - v0.2.0 - Require Julia v1.9 for extensions - Create package extensions diff --git a/Project.toml b/Project.toml index fd944b59..48abe72c 100644 --- a/Project.toml +++ b/Project.toml @@ -1,7 +1,7 @@ name = "EcoSISTEM" uuid = "ed2dc23b-ada4-5fdb-a26f-56368a14ad8f" authors = ["Claire Harris ", "Richard Reeve "] -version = "0.2.0" +version = "0.2.1" [deps] ArchGDAL = "c9ce4bd3-c3d5-55b8-8973-c0e20141b8c3" diff --git a/src/EcoSISTEM.jl b/src/EcoSISTEM.jl index 9a74ea26..70ebc862 100644 --- a/src/EcoSISTEM.jl +++ b/src/EcoSISTEM.jl @@ -95,16 +95,18 @@ export meta_simpson, meta_shannon, meta_speciesrichness, mean_abun, using Random abstract type MPIGridLandscape end +export MPIGridLandscape + abstract type MPIEcosystem{MPIGL <: MPIGridLandscape, Part <: AbstractAbiotic, SL <: SpeciesList, TR <: AbstractTraitRelationship} <: AbstractEcosystem{Part, SL, TR} end - -export MPIGridLandscape export MPIEcosystem + function gather_abundance end function gather_diversity end export gather_abundance, gather_diversity + function emptyMPIgridlandscape end function synchronise_from_rows! end function synchronise_from_cols! end