template class BRKGA::BRKGA_MP_IPR

Overview

This class represents a Multi-Parent Biased Random-key Genetic Algorithm with Implicit Path Relinking (BRKGA-MP-IPR). More…

#include <brkga_mp_ipr.hpp>

template <class Decoder>
class BRKGA_MP_IPR {
public:
    // construction

    BRKGA_MP_IPR(
        Decoder& decoder_reference,
        const Sense sense,
        const unsigned seed,
        const unsigned chromosome_size,
        const BrkgaParams& params,
        const unsigned max_threads = 1,
        const bool evolutionary_mechanism_on = true
    );

    // methods

    void setInitialPopulation(const std::vector<Chromosome>& chromosomes);
    void setBiasCustomFunction(const std::function<double(const unsigned)>& func);
    void setShakingMethod(const std::function<void(double lower_bound, double upper_bound, std::vector<std::shared_ptr<Population>>&populations, std::vector<std::pair<unsigned, unsigned>>&shaken)>& func);
    void setStoppingCriteria(const std::function<bool(const AlgorithmStatus&)>& stopping_criteria);
    void addNewSolutionObserver(const std::function<bool(const AlgorithmStatus&)>& func);

    AlgorithmStatus run(
        const ControlParams& control_params,
        std::ostream* logger = &std::cout
    );

    void evolve(unsigned generations = 1);

    PathRelinking::PathRelinkingResult pathRelink(
        PathRelinking::Type pr_type,
        PathRelinking::Selection pr_selection,
        std::shared_ptr<DistanceFunctionBase> dist,
        unsigned number_pairs,
        double minimum_distance,
        std::size_t block_size = 1,
        std::chrono::seconds max_time = std::chrono::seconds{0},
        double percentage = 1.0
    );

    PathRelinking::PathRelinkingResult pathRelink(
        std::shared_ptr<DistanceFunctionBase> dist,
        std::chrono::seconds max_time = std::chrono::seconds{0}
    );

    void exchangeElite(unsigned num_immigrants);
    void reset();

    void shake(
        unsigned intensity,
        ShakingType shaking_type,
        unsigned population_index = std::numeric_limits<unsigned>::infinity()
    );

    void injectChromosome(
        const Chromosome& chromosome,
        unsigned population_index,
        unsigned position
    );

    const Population& getCurrentPopulation(unsigned population_index = 0) const;
    const Chromosome& getBestChromosome() const;
    fitness_t getBestFitness() const;

    const Chromosome& getChromosome(
        unsigned population_index,
        unsigned position
    ) const;

    fitness_t getFitness(
        unsigned population_index,
        unsigned position
    ) const;

    const BrkgaParams& getBrkgaParams() const;
    Sense getOptimizationSense() const;
    unsigned getChromosomeSize() const;
    unsigned getEliteSize() const;
    unsigned getNumMutants() const;
    bool evolutionaryIsMechanismOn() const;
    unsigned getMaxThreads() const;
};

Detailed Documentation

This class represents a Multi-Parent Biased Random-key Genetic Algorithm with Implicit Path Relinking (BRKGA-MP-IPR).

Carlos Eduardo de Andrade ce.andrade@gmail.com

2023

Main capabilities

Evolutionary process

In the BRKGA-MP-IPR, we keep a population of chromosomes divided between the elite and the non-elite group. During the mating, multiple parents are chosen from the elite group and the non-elite group. They are sorted either on no-decreasing order for minimization or non-increasing order to maximization problems. Given this order, a bias function is applied to the rank of each chromosome, resulting in weight for each one. Using a roulette method based on the weights, the chromosomes are combined using a biased crossover.

This code also implements the island model, where multiple populations can be evolved in parallel, and migration between individuals between the islands are performed using exchangeElite() method.

This code requires the template argument Decoder be a class or functor object capable to map a chromosome to a solution for the specific problem, and return a value to be used as fitness to the decoded chromosome. The decoder must have the method

double decode(Chromosome& chromosome, bool rewrite);

where Chromosome is a vector<double> (in general) representing a solution and rewrite is a boolean indicating that if the decode should rewrite the chromosome in case it implements local searches and modifies the initial solution decoded from the chromosome. Since this API has the capability of decoding several chromosomes in parallel, the user must guarantee that Decoder::decode(...) is thread-safe. Therefore, we do recommend to have the writable variables per thread. Please, see the example that follows this code.

Implicit Path Relinking

This API also implements the Implicit Path Relinking leveraging the decoder capabilities. To perform the path relinking, the user must call pathRelink() method, indicating the type of path relinking (direct or permutation-based, see PathRelinking::Type), the selection criteria (best solution or random elite, see PathRelinking::Selection), the distance function (to choose individuals far enough, see BRKGA::DistanceFunctionBase, BRKGA::HammingDistance, and BRKGA::KendallTauDistance), a maximum time or a maximum path size.

In the presence of multiple populations, the path relinking is performed between elite chromosomes from different populations, in a circular fashion. For example, suppose we have 3 populations. The framework performs 3 path relinkings: the first between individuals from populations 1 and 2, the second between populations 2 and 3, and the third between populations 3 and 1. In the case of just one population, both base and guiding individuals are sampled from the elite set of that population.

Note that the algorithm tries to find a pair of base and guiding solutions with a minimum distance given by the distance function. If this is not possible, a new pair of solutions are sampled (without replacement) and tested against the distance. In case it is not possible to find such pairs for the given populations, the algorithm skips to the next pair of populations (in a circular fashion, as described above). Yet, if such pairs are not found in any case, the algorithm declares failure. This indicates that the populations are very homogeneous.

The API will call Decoder::decode() always with rewrite = false. The reason is that if the decoder rewrites the chromosome, the path between solutions is lost and inadvertent results may come up. Note that at the end of the path relinking, the method calls the decoder with rewrite = true in the best chromosome found to guarantee that this chromosome is re-written to reflect the best solution found.

Other capabilities

Multi-start

This API also can be used as a simple multi-start algorithm without evolution. To do that, the user must supply in the constructor the argument evolutionary_mechanism_on = false. This argument makes the elite set has one individual and the number of mutants n - 1, where n is the size of the population. This setup disables the evolutionary process completely.

Initial Population

This API allows the user provides a set of initial solutions to warm start the algorithm. In general, initial solutions are created using other (fast) heuristics and help the convergence of the BRKGA. To do that, the user must encode the solutions on Chromosome (vector<double>) and pass to the method setInitialPopulation() as a vector<#Chromosome>.

General Usage

  1. The user must call the BRKGA_MP_IPR constructor and pass the desired parameters. Please, see BRKGA_MP_IPR::BRKGA_MP_IPR for parameter details;

  2. (Optional) The user provides the warm start solutions using setInitialPopulation();

  3. (Optional) The user provides a callback to track the optimization using addNewSolutionObserver();

  4. (Optional) The user provides custom stopping criteria function using setStoppingCriteria();

  5. (Optional) The user provides a custom shaking procedure using setShakingMethod();

  6. The user calls run() to cary out the optimization; This method returns a AlgorithmStatus object with the results of the optimization.

For a comprehensive and detailed usage, please see the examples that follow this API.

About multi-threading

This API is capable of decoding several chromosomes in parallel, as mentioned before. This capability is based on OpenMP (http://www.openmp.org) and the compiler must have support to it. Most recent versions of GNU G++ and Intel C++ Compiler support OpenMP. Clang supports OpenMP since 3.8. However, there are some issues with the libraries, and sometimes, the parallel sections are not enabled. On the major, the user can find fixes to his/her system.

Since, in general, the decoding process can be complex and lengthy, it is recommended that the number of threads used DO NOT exceed the number of physical cores in the machine. This improves the overall performance drastically, avoiding cache misses and racing conditions. Note that the number of threads is also tied to the memory utilization, and it should be monitored carefully.

History

This API was based on the code by Rodrigo Franco Toso, Sep 15, 2011. http://github.com/rfrancotoso/brkgaAPI

Construction

BRKGA_MP_IPR(
    Decoder& decoder_reference,
    const Sense sense,
    const unsigned seed,
    const unsigned chromosome_size,
    const BrkgaParams& params,
    const unsigned max_threads = 1,
    const bool evolutionary_mechanism_on = true
)

Builds the algorithm and its data strtuctures with the given arguments.

Parameters:

decoder_reference

a reference to the decoder object. NOTE: BRKGA uses such object directly for decoding.

sense

the optimization sense (maximization or minimization).

seed

the seed for the random number generator.

chromosome_size

number of genes in each chromosome.

params

BRKGA and IPR parameters object loaded from a configuration file or manually created. All the data is copied.

max_threads

number of threads to perform parallel decoding.

NOTE : Decoder::decode() MUST be thread-safe.

evolutionary_mechanism_on

if false, no evolution is performed but only chromosome decoding. Very useful to emulate a multi-start algorithm.

std::range_error

if some parameter or combination of parameters does not fit.

Methods

void setInitialPopulation(const std::vector<Chromosome>& chromosomes)

Sets individuals to initial population.

Set initial individuals into the population to work as warm-starters. Such individuals can be obtained from solutions of external procedures such as fast heuristics, other methaheuristics, or even relaxations from a mixed integer programming model that models the problem.

We assign as many individuals as possible across all populations. Extra individuals are disregarded.

Parameters:

chromosomes

a set of individuals encoded as Chromosomes.

void setBiasCustomFunction(const std::function<double(const unsigned)>& func)

Sets a custom bias function used to build the probabilities.

It must be a positive non-increasing function, i.e. \(f: \mathbb{N}^+ \to \mathbb{R}^+\) such that \(f(i) \ge 0\) and \(f(i) \ge f(i+1)\) for \(i \in [1, \ldots, total\_parents]\). For example

setBiasCustomFunction(
    [](const unsigned x) {
        return 1.0 / (x * x);
    }
);

sets an inverse quadratic function.

Parameters:

func

a reference to a unary positive non-increasing function.

std::runtime_error

in case the function is not a non-decreasing positive function.

void setShakingMethod(const std::function<void(double lower_bound, double upper_bound, std::vector<std::shared_ptr<Population>>&populations, std::vector<std::pair<unsigned, unsigned>>&shaken)>& func)

Sets a custom shaking procedure.

For more details, see BrkgaParams::custom_shaking.

Parameters:

func

a callback function. For example, the code below implements the standard mutation:

// A random number generator.
std::mt19937 rng(2700001);
rng.discard(rng.state_size);

// Change some values from elite chromosomes from all populations.
// Similar to a standard mutation.
algorithm.setShakingMethod(
    [&](double lower_bound, double upper_bound,
        std::vector<std::shared_ptr<Population>>& populations,
        std::vector<std::pair<unsigned, unsigned>>& shaken) {

        // Determines whether we change the allele or not.
        std::bernoulli_distribution must_change(0.50);

        // Determines the value of the allele.
        std::uniform_real_distribution<> allele_value(lower_bound, upper_bound);

        for(unsigned pop_idx = 0; pop_idx < populations.size(); ++pop_idx) {
            auto& population = populations[0]->population;
            for(unsigned chr_idx = 0; chr_idx < population.size(); ++chr_idx) {
                auto& chromosome = population[chr_idx];

                bool change = false;
                for(unsigned i = 0; i < chromosome.size(); ++i) {
                    if(must_change(rng)) {
                        chromosome[i] = allele_value(rng);
                        change = true;
                    }
                }

                if(change)
                    shaken.push_back({pop_idx, chr_idx});
            }
        }
    };
);
void setStoppingCriteria(const std::function<bool(const AlgorithmStatus&)>& stopping_criteria)

Sets a custom stopping criteria supplied by the user.

The algorithm always test for the maximum running time and for the maximum stalled iterations/generations given by ControlParams indenpendently of the stopping criteria function supplied by the user. This is especially important when activating the implicit path reliking which is very timing consuming.

Warning

If you are using IPR, we STRONGLY RECOMMEND TO SET A MAXIMUM TIME since this is the core stopping criteria on IPR.

If you really mean to have no maximum time and/or maximum stalled iterations set, we recommend to use the following code:

// After reading your parameters, e.g.,
// auto [brkga_params, control_params] = readConfiguration("config.conf");

// You can set to the max.
control_params.maximum_running_time = std::chrono::seconds::max();

control_params.stall_offset = numeric_limits<unsigned>::max();

Parameters:

stopping_criteria

a callback function to determine is the algorithm must stop. For instance, the following lambda function tests if the best solution reached a given value for a minimization problem:

fitness_t my_magical_solution = 10;

algorithm.setStoppingCriteria(
    [&](const AlgorithmStatus& status) {
        return status.best_fitness <= my_magical_solution;
    }
);
void addNewSolutionObserver(const std::function<bool(const AlgorithmStatus&)>& func)

Adds a callback function called when the best solution is improved.

It must take a reference to AlgorithmStatus and return true if the algorithm should stop immediately. You may have as many observers as you want. They will be called in the order they are added.

Parameters:

func

a callback function such as

bool check_solution(const AlgorithmStatus& status) {
    std::cout << "\n" << status.best_fitness;
    return true; // Stop the optimization.
}
//...
algorithm.addNewSolutionObserver(check_solution);

or a lambda function such as

algorithm.addNewSolutionObserver(
    [](const AlgorithmStatus& status) {
        std::cout
        << "> Iter: " << status.current_iteration
        << " | solution: " << status.best_fitness
        << " | time: " << status.current_time
        << std::endl;
        return false; // Dont' stop the optimization.
     }
);
AlgorithmStatus run(
    const ControlParams& control_params,
    std::ostream* logger = &std::cout
)

Runs the full framework performing evolution, path-reliking, exchanges, shakes, and resets according to the parameters.

This method uses all facilities associated with this BRKGA-MP-IPR library, providing a comprehensive and easy-to-use single-entry point. The main loop always evolves one generation per iteration and calls other procedures based on the number of stalled iterations (i.e., the number of iterations without improvement in the best solution), and the given user thresholds in Control Parameters. If the thresholds are all the same, the main loop should be like this:

while(!must_stop) {
    evolve(); // One generation.
    if(best solution improvement) {
        Save best solution;
        Call observer callbacks;
    }

    if(!must_stop && ipr_interval > 0 && stalled_iterations > 0 &&
       stalled_iterations % ipr_interval == 0) {
         pathRelink();
         if(best solution improvement) {
             Save best solution;
             Call observer callbacks;
         }
     }

    if(!must_stop && exchange_interval > 0 && stalled_iterations > 0 &&
       stalled_iterations % exchange_interval == 0) {
         exchangeElite();
    }

    if(!must_stop && shake_interval > 0 && stalled_iterations > 0 &&
       stalled_iterations % shake_interval == 0) {
         shake();
    }

    if(!must_stop && reset_interval > 0 && stalled_iterations > 0 &&
       stalled_iterations % reset_interval == 0) {
         reset();
    }
}

Therefore, note that the order that pathRelink(), exchangeElite(), shake(), and reset() are called, depends on the thresholds defined in ControlParams.

For path relinking, the block size is computed by \(\lceil \alpha \times \sqrt{p} \rceil\) where \(\alpha\) is BrkgaParams::alpha_block_size and \(p\) is BrkgaParams::population_size. If the size is larger than the chromosome size, the size is set to half of the chromosome size. For more details on that, refer to pathRelink().

Note

The algorithm always test against maximum running time and for the maximum stalled iterations/generations given by ControlParams indenpendently of the stopping criteria function supplied by the user. This is especially important when activating the implicit path reliking which is very timing consuming. If you are using IPR, we STRONGLY RECOMMEND TO SET A MAXIMUM TIME since this is the core stopping criteria on IPR.

Warning

The decoding is done in parallel using threads, and the user must guarantee that the decoder is THREAD-SAFE. If such property cannot be held, we suggest using a single thread for optimization.

Parameters:

control_params

the parameters to control the algorithm flow, such as calling exchanges, shakes, and IPR.

logger

a output stream to log some information.

std::runtime_error

in the following cases:

  1. IPR is active (ipr_interva > 0) but the distance function is not set;

  2. Shaking is active (shake_interval > 0) and it is set as ‘CUSTOM’. However the custom shaking procedure was not supplied.

  3. Shaking is active (shake_interval > 0). However, the intensity bounds are out of range. Should be (0.0, 1.0] and ‘shaking_intensity_lower_bound <= shaking_intensity_upper_bound’.

Returns:

The last algorithm status before the stopping criteria are met.

void evolve(unsigned generations = 1)

Evolves the current populations following the guidelines of Multi-parent BRKGAs.

Warning

The decoding is done in parallel using threads, and the user must guarantee that the decoder is THREAD-SAFE. If such property cannot be held, we suggest using a single thread for optimization.

Parameters:

generations

number of generations to be evolved. Must be larger than zero.

std::range_error

if the number of generations is zero.

PathRelinking::PathRelinkingResult pathRelink(
    PathRelinking::Type pr_type,
    PathRelinking::Selection pr_selection,
    std::shared_ptr<DistanceFunctionBase> dist,
    unsigned number_pairs,
    double minimum_distance,
    std::size_t block_size = 1,
    std::chrono::seconds max_time = std::chrono::seconds{0},
    double percentage = 1.0
)

Performs path relinking between elite solutions that are, at least, a given minimum distance between themselves. In this method, the local/loaded parameters are ignored in favor to the supplied ones.

In the presence of multiple populations, the path relinking is performed between elite chromosomes from different populations, in a circular fashion. For example, suppose we have 3 populations. The framework performs 3 path relinkings: the first between individuals from populations 1 and 2, the second between populations 2 and 3, and the third between populations 3 and 1. In the case of just one population, both base and guiding individuals are sampled from the elite set of that population.

Note that the algorithm tries to find a pair of base and guiding solutions with a minimum distance given by the distance function. If this is not possible, a new pair of solutions are sampled (without replacement) and tested against the distance. In case it is not possible to find such pairs for the given populations, the algorithm skips to the next pair of populations (in a circular fashion, as described above). Yet, if such pairs are not found in any case, the algorithm declares failure. This indicates that the populations are very homogeneous.

If the found solution is the best solution found so far, IPR replaces the worst solution by it. Otherwise, IPR computes the distance between the found solution and all other solutions in the elite set, and replaces the worst solution by it if and only if the found solution is, at least, minimum_distance from all them.

The API will call Decoder::decode() function always with rewrite = false. The reason is that if the decoder rewrites the chromosome, the path between solutions is lost and inadvertent results may come up. Note that at the end of the path relinking, the method calls the decoder with rewrite = true in the best chromosome found to guarantee that this chromosome is re-written to reflect the best solution found.

This method is a multi-thread implementation. Instead of to build and decode each chromosome one at a time, the method builds a list of candidates, altering the alleles/keys according to the guide solution, and then decode all candidates in parallel. Note that O(chromosome_size^2 / block_size) additional memory is necessary to build the candidates, which can be costly if the chromosome_size is very large.

Warning

As it is in evolve(), the decoding is done in parallel using threads, and the user must guarantee that the decoder is THREAD-SAFE. If such property cannot be held, we suggest using a single thread for optimization.

Parameters:

pr_type

type of path relinking to be performed. See PathRelinking::Type.

pr_selection

of which individuals use to path relinking. See PathRelinking::Selection.

dist

a pointer to a functor/object to compute the distance between two chromosomes. This object must be inherited from BRKGA::DistanceFunctionBase and implement its methods.

number_pairs

number of chromosome pairs to be tested. If 0, all pairs are tested.

minimum_distance

between two chromosomes computed by dist.

block_size

number of alleles to be exchanged at once in each iteration. If one, the traditional path relinking is performed.

max_time

aborts path relinking when reach max_time. If max_time <= 0, no limit is imposed.

percentage

defines the size, in percentage, of the path to build. Default: 1.0 (100%).

std::range_error

if the percentage or size of the path is not in (0, 1].

Returns:

A PathRelinking::PathRelinkingResult depending on the relink status.

PathRelinking::PathRelinkingResult pathRelink(
    std::shared_ptr<DistanceFunctionBase> dist,
    std::chrono::seconds max_time = std::chrono::seconds{0}
)

Performs path relinking between elite solutions that are, at least, a given minimum distance between themselves.

This method uses all parameters supplied in the constructor. In particular, the block size is computed by \(\lceil \alpha \times \sqrt{p} \rceil\) where \(\alpha\) is BrkgaParams::alpha_block_size and \(p\) is BrkgaParams::population_size. If the size is larger than the chromosome size, the size is set to half of the chromosome size.

Please, refer to pathRelink() for details.

Parameters:

dist

a pointer to a functor/object to compute the distance between two chromosomes. This object must be inherited from BRKGA::DistanceFunctionBase and implement its methods.

max_time

aborts path relinking when reach max_time. If max_time <= 0, no limit is imposed.

std::range_error

if the percentage or size of the path is not in (0, 1].

Returns:

A PathRelinking::PathRelinkingResult depending on the relink status.

void exchangeElite(unsigned num_immigrants)

Exchanges elite-solutions between the populations.

Given a population, the num_immigrants best solutions are copied to the neighbor populations, replacing their worth solutions. If there is only one population, nothing is done.

Parameters:

num_immigrants

number of elite chromosomes to select from each population.

std::range_error

if the number of immigrants less than one or it is larger than or equal to the population size divided by the number of populations minus one, i.e. \(\lceil \frac{population\_size}{num\_independent\_populations} \rceil - 1\).

void reset()

Resets all populations with brand new keys.

All warm-start solutions provided setInitialPopulation() are discarded. You may use injectChromosome() to insert those solutions again.

void shake(
    unsigned intensity,
    ShakingType shaking_type,
    unsigned population_index = std::numeric_limits<unsigned>::infinity()
)

Performs a shaking in the chosen population.

Parameters:

intensity

the intensity of the shaking.

shaking_type

either CHANGE or SWAP moves.

population_index

the index of the population to be shaken. If population_index >= num_independent_populations, all populations are shaken.

void injectChromosome(
    const Chromosome& chromosome,
    unsigned population_index,
    unsigned position
)

Injects/Replaces a chromosome of a given population into a given position.

The new chromosome replaces the old one, and the decoder is triggered to compute the new fitness. Once done, the population is re-sorted according to the chromosomes’ fitness.

Parameters:

chromosome

the chromosome to be injected.

population_index

the population index.

position

the chromosome position.

std::range_error

eitheir if population_index is larger than number of populations; or position is larger than the population size; or chromosome.size() != chromosome_size

const Population& getCurrentPopulation(unsigned population_index = 0) const

Returns a reference to a current population.

Parameters:

population_index

the population index.

std::range_error

if the index is larger than number of populations.

const Chromosome& getBestChromosome() const

Returns a reference to the chromosome with best fitness among all current populations.

Warning

Note that this method does not return the best solution but the one within the current population. If a shake() or reset() is called, the best solution may be lose in the populations. However, if you are using run(), the best solution is returned by that method. If not, you must keep track of the best solution.

fitness_t getBestFitness() const

Returns the best fitness among all current populations.

Warning

Note that this method does not return the best fitness but the one within the current population. If a shake() or reset() is called, the best fitness may be lose in the populations. However, if you are using run(), the best fitness is returned by that method. If not, you must keep track of the best fitness.

const Chromosome& getChromosome(
    unsigned population_index,
    unsigned position
) const

Returns a reference to a chromosome of the given population.

Parameters:

population_index

the population index.

position

the chromosome position, ordered by fitness. The best chromosome is located in position 0.

std::range_error

eitheir if population_index is larger than number of populations, or position is larger than the population size.

fitness_t getFitness(
    unsigned population_index,
    unsigned position
) const

Returns the fitness of a chromosome of the given population.

Parameters:

population_index

the population index.

position

the chromosome position, ordered by fitness. The best chromosome is located in position 0.

std::range_error

eitheir if population_index is larger than number of populations, or position is larger than the population size.