GPLIB++
mtanisoga.cpp
Go to the documentation of this file.
1 #include <iostream>
2 #include <fstream>
3 #include <algorithm>
4 #include <numeric>
5 #include <sstream>
6 #include <string>
7 #include "SimpleSelect.h"
9 #include "UniformRNG.h"
10 #include "gentypes.h"
11 #include "BinaryPopulation.h"
12 #include "StandardPropagation.h"
13 #include "GrayTranscribe.h"
14 #include "AnnealingGA.h"
15 #include "ParetoGA.h"
16 #include "Aniso1DMTObjective.h"
17 #include "PTensor1DMTObjective.h"
18 #include "MTStation.h"
19 #include "Adaptors.h"
20 #include "MTFitSetup.h"
21 #include "PTensorMTStation.h"
22 #include <boost/bind.hpp>
23 #include <boost/shared_ptr.hpp>
24 #include <boost/filesystem.hpp>
25 #include <boost/date_time/posix_time/posix_time.hpp>
26 #include <complex>
27 #include "Util.h"
28 #include "MTAnisoRoughness.h"
29 #include "MTAnisoGAConf.h"
30 
31 using namespace std;
32 using namespace gplib;
33 
34 enum tgatype
35  {
37  };
38 typedef boost::shared_ptr<GeneralObjective> pCGeneralObjective;
39 
40 
41 
42 void SetupAnnealingGA(boost::shared_ptr<GeneralGA> &GA,
43  const MTAnisoGAConf Configuration)
44  {
45  double InitTemperature = Configuration.inittemp;
46  if (Configuration.inittemp < 0)
47  {
48  double avgcombfit = 0;
49  GA->DoIteration(0, false);
50  /*for (int i = 0; i < GA->GetAvgFit().size(); ++i)
51  avgcombfit += GA->GetAvgFit().at(i);*/
52  avgcombfit = accumulate(GA->GetAvgFit().begin(), GA->GetAvgFit().end(),
53  0.0);
54  InitTemperature = avgcombfit * abs(Configuration.inittemp);
55  }
56  AnnealingGA *AnnealGA = dynamic_cast<AnnealingGA*> (GA.get());
57  AnnealGA->SetParams(InitTemperature, Configuration.annealinggeneration,
58  Configuration.coolingratio);
59  }
60 
61 //! Program to invert MT data for 1D anisotropic structure with a genetic algorithm
62 int main(int argc, char *argv[])
63  {
64  boost::posix_time::ptime starttime =
65  boost::posix_time::microsec_clock::local_time();
66  string version = "$Id: mtanisoga.cpp 1849 2010-05-07 11:53:45Z mmoorkamp $";
67  cout << endl << endl;
68  cout << "Program " << version << endl;
69  cout << "Genetic algorithm to jointly invert MT data for 1D anisotropy"
70  << endl;
71  cout << "look at mtanisoga.conf for configuration and setup " << endl;
72  cout << "The root name of log-files can be overwritten by using " << endl;
73  cout << "1dinvga newrootname " << endl;
74 
75  MTAnisoGAConf Configuration;
76 
77  UniformRNG Random;
78  try
79  {
80  Configuration.GetData("mtanisoga.conf");
81 
82  ttranscribed transcribed;
83  tgatype GAtype = pareto;
84 
85  string outputbase;
86 
87  if (argc > 1)
88  outputbase = argv[1];
89  else
90  outputbase = Configuration.outputbase;
91 
92  ofstream misfitfile;
93  ofstream valuefile;
94  ofstream probfile;
95  ofstream distancefile;
96  ofstream fitstat((outputbase + ".ftst").c_str());
97  ofstream uniquepop((outputbase + ".unpop").c_str());
98  ofstream rankfile((outputbase + ".rank").c_str());
99  ofstream frontfile((outputbase + ".front").c_str());
100  ofstream bestfile((outputbase + ".best").c_str());
101 
102  if (Configuration.verbose) // in most cases we do not need all of these values
103  {
104  misfitfile.open((outputbase + ".msft").c_str());
105  valuefile.open((outputbase + ".vals").c_str());
106  probfile.open((outputbase + ".probs").c_str());
107  distancefile.open((outputbase + ".dist").c_str());
108  }
109 
110  if ((Configuration.thickbase.size() != Configuration.resbase.size())
111  || (Configuration.thickbase.size()
112  != Configuration.anisobase.size())
113  || (Configuration.thickbase.size()
114  != Configuration.strikebase.size()))
115  {
116  cerr << "Configuration of genes is not consistent ! ";
117  return 100;
118  }
119  const int nparam = 4;
120  const int genes = 4 * Configuration.thickbase.size();
121  const int popsize = Configuration.popsize;
122  const int maxgen = Configuration.generations;
123  const int paramlength = genes / nparam;
124  ttranscribed basevalues(genes);
125  ttranscribed stepsizes(genes);
126  tsizev genesizes(genes);
127 
128  misfitfile << popsize << " " << maxgen << " " << endl;
129  valuefile << genes << " " << popsize << " " << maxgen << endl;
130 
131  copy(Configuration.resbase.begin(), Configuration.resbase.end(),
132  basevalues.begin());
133  copy(Configuration.resstep.begin(), Configuration.resstep.end(),
134  stepsizes.begin());
135  copy(Configuration.ressizes.begin(), Configuration.ressizes.end(),
136  genesizes.begin());
137 
138  copy(Configuration.thickbase.begin(), Configuration.thickbase.end(),
139  basevalues.begin() + paramlength);
140  copy(Configuration.thickstep.begin(), Configuration.thickstep.end(),
141  stepsizes.begin() + paramlength);
142  copy(Configuration.thicksizes.begin(), Configuration.thicksizes.end(),
143  genesizes.begin() + paramlength);
144 
145  copy(Configuration.anisobase.begin(), Configuration.anisobase.end(),
146  basevalues.begin() + 2 * paramlength);
147  copy(Configuration.anisostep.begin(), Configuration.anisostep.end(),
148  stepsizes.begin() + 2 * paramlength);
149  copy(Configuration.anisosizes.begin(), Configuration.anisosizes.end(),
150  genesizes.begin() + 2 * paramlength);
151 
152  copy(Configuration.strikebase.begin(), Configuration.strikebase.end(),
153  basevalues.begin() + 3 * paramlength);
154  copy(Configuration.strikestep.begin(), Configuration.strikestep.end(),
155  stepsizes.begin() + 3 * paramlength);
156  copy(Configuration.strikesizes.begin(),
157  Configuration.strikesizes.end(), genesizes.begin() + 3
158  * paramlength);
159 
160  const int totalgenesize = accumulate(genesizes.begin(),
161  genesizes.end(), 0);
162  BinaryPopulation Population(popsize, totalgenesize, Random, true);
163  BinaryTournamentSelect Select(Random, boost::bind(
164  &GeneralPopulation::GetProbabilities, boost::ref(Population)),
165  boost::bind(&GeneralPopulation::GetCrowdingDistances, boost::ref(
166  Population)));
167  StandardPropagation Propagator(&Select, &Population, &Random);
168  double indmutationprob = 1.0 - std::pow(1.0 - Configuration.mutationprob, 1.
169  / totalgenesize); // we specify the overall probability for mutation in config file
170  cout << "Individual Mutation set to: " << indmutationprob << endl;
171  Propagator.SetParams(indmutationprob, Configuration.crossoverprob);
172  GrayTranscribe Transcribe(basevalues, stepsizes, genesizes);
173 
174  //read in data
175 
176 
177  //setup MT Objective function
178  boost::shared_ptr<PlottableObjective> MTObjective;
179 
180  if (Configuration.mtfit == "ptensor")
181  {
182  PTensorMTStation PTData;
183  PTData.GetData(Configuration.ptensordata);
184  boost::shared_ptr<PTensor1DMTObjective> PTObjective(
185  new PTensor1DMTObjective(PTData));
186  PTObjective->SetErrorLevel(Configuration.phaseerror);
187  MTObjective = PTObjective;
188  }
189  else
190  {
192  MTData.GetData(Configuration.mtinputdata);
193  boost::shared_ptr<Aniso1DMTObjective> AnisoObjective(
194  new Aniso1DMTObjective(MTData));
195  SetupMTFitParameters(Configuration, *AnisoObjective);
196  MTObjective = AnisoObjective;
197  }
198  MTObjective->SetFitExponent(Configuration.mtfitexponent);
199 
200  //setup regularization
201 
202  boost::shared_ptr<MTAnisoRoughness> MTRegularization(
203  new MTAnisoRoughness());
204  MTRegularization->SetCondDiffWeight(Configuration.conddiffweight);
205  MTRegularization->SetAnisotropyWeight(Configuration.anisotropyweight);
206  MTRegularization->SetStrikeDiffWeight(Configuration.strikediffweight);
207  vector<pCGeneralObjective> ObjVector;
208  ObjVector.push_back(MTObjective);
209  ObjVector.push_back(MTRegularization);
210 
211  boost::shared_ptr<GeneralGA> GA;
212  if (Configuration.gatype == "pareto")
213  {
214  GAtype = pareto;
215  boost::shared_ptr<ParetoGA> Pareto(new ParetoGA(&Propagator,
216  &Population, &Transcribe, ObjVector, Configuration.threads));
217  GA = Pareto;
218  }
219  else if (Configuration.gatype == "anneal")
220  {
221  GAtype = anneal;
222  boost::shared_ptr<AnnealingGA> Annealing(new AnnealingGA(
223  &Propagator, &Population, &Transcribe, ObjVector,
224  Configuration.threads));
225  GA = Annealing;
226  }
227  else
228  {
229  cerr << "Invalid GA type in configuration file." << endl;
230  cerr << " Has to be either 'pareto' or 'anneal' !" << endl;
231  return 100;
232  }
233  GA->SetWeights(Configuration.weights);
234  //Setting up indices for forward modelling
235  GeneralGA::tparamindv ParamIndices;
236  const int nmtparam = 4 * Configuration.resbase.size();
237  std::vector<int> LocalIndices;
238  generate_n(back_inserter(LocalIndices), nmtparam, IntSequence(0)); //generate Indices 0 .. nmtparam for MT Objective function
239  ParamIndices.push_back(LocalIndices);//both data misfit and regularization
240  ParamIndices.push_back(LocalIndices); //work on all parts of the model vector
241 
242  GA->SetParameterIndices(ParamIndices);
243  GA->SetElitist(Configuration.elitist);
244  //we have to do a few things that are specific to the type of genetic algorithm
245  if (GAtype == anneal)
246  {
247  SetupAnnealingGA(GA, Configuration);
248  Population.InitPop();
249  }
250  //now comes the core loop
251  for (int i = 0; i < maxgen; ++i) //iterate over the specified number of generattions
252  {
253  GA->DoIteration(i, i == (maxgen - 1)); // pass iteration number to GA and tell whether it's the final iteration
254  if (Configuration.verbose) // in most cases we do not need all of these values
255  {
256  GA->PrintMisfit(misfitfile);
257  GA->PrintTranscribed(valuefile);
258  Population.PrintProbabilities(probfile);
259  Population.PrintDistances(distancefile);
260  }
261  GA->PrintBestMisfit(frontfile);
262 
263  fitstat << i << " ";
264  GA->PrintFitStat(fitstat);
265  }
266 
267  //the following is to organize the output of the best models
268  int bestcount = 0;
269  UniquePop FinalUnique; //we only want to write each model once, even if it is several times in the population
270  tfitvec bestfit(ObjVector.size());
271  tpopmember member(Population.GetGenesize());
272 
273  vector<int> BestIndices(GA->GetBestModelIndices()); // store the indices of the best models
274  const unsigned int noutmodels = GA->GetNBestmodels(); // and how many of them we have
275  assert(BestIndices.size() == noutmodels);
276  const unsigned int nobjective = ObjVector.size();
277  GeneralGA::tparamvector LocalParameters(nobjective);
278  for (unsigned int i = 0; i < nobjective; ++i)
279  LocalParameters.at(i).resize(ParamIndices.at(i).size(), false);
280  ofstream finalmodels((outputbase + ".final").c_str());
281  for (unsigned int h = 0; h < noutmodels; ++h) //for each best model
282  {
283  member = row(Population.GetPopulation(), BestIndices.at(h)); //get the genetic string
284 
285 
286  transcribed = Transcribe.GetValues(member); //transcribe the genetic string to model parameters
287  GA->SetupParams(transcribed, LocalParameters); //copy the right parameters for each objective function
288 
289  for (unsigned int f = 0; f < nobjective; ++f) //calculate the performance for each function
290  {
291  if (Configuration.weights.at(f) > 0)
292  bestfit(f) = ObjVector.at(f)->CalcPerformance(
293  LocalParameters.at(f));
294  else
295  bestfit(f) = 0;
296 
297  }
298 
299  if (FinalUnique.Insert(bestfit, transcribed)) //if we can insert it, it is a new model that we didn't output yet
300  {
301  copy(bestfit.begin(), bestfit.end(), ostream_iterator<double> (
302  finalmodels, " ")); //print the misfits in the .final file
303  finalmodels << " " << bestcount << endl; // write the index in the last column for easier identification
304  ostringstream filename;
305  filename << "best_" << bestcount;
306  cout << "Best: " << h << " Index : " << BestIndices.at(h)
307  << " "; //write some info to screen
308  for (unsigned int f = 0; f < nobjective; ++f)
309  cout << bestfit(f) << " ";
310  cout << " Saved as : " << filename.str();
311 
312  MTObjective->WriteModel(filename.str() + "_mt.mod"); //write model, plotting file and synthetic data
313  MTObjective->WritePlot(filename.str() + "_mt.plot");
314  MTObjective->WriteData(filename.str());
315 
316  filename.str("");
317  bestcount++;
318  cout << endl;
319  }
320 
321  }
322 
323  if (Configuration.verbose) // in most cases we do not need all of these values
324  {
325  GA->PrintUniquePop(uniquepop);
326  FinalUnique.PrintAll(bestfile);
327  }
328  boost::posix_time::ptime endtime =
329  boost::posix_time::microsec_clock::local_time();
330  std::cout << "Total Runtime: " << endtime - starttime << std::endl;
331  } catch (FatalException &e)
332  {
333  cerr << e.what() << endl;
334  return -1;
335  }
336  }
void SetupMTFitParameters(const tConfObject &Configuration, C1DMTObjective &Objective)
Definition: MTFitSetup.h:8
ublas::vector< double > ttranscribed
Definition: gentypes.h:21
Generate uniformly distributed random numbers, this is basically a wrapper for the boost random numbe...
Definition: UniformRNG.h:19
bool Insert(const tfitvec &fitness, const ttranscribed &popmember)
Definition: UniquePop.cpp:17
Caclulate the roughness for anisotropic MT models.
boost::shared_ptr< GeneralObjective > pCGeneralObjective
Definition: mtanisoga.cpp:38
void SetupAnnealingGA(boost::shared_ptr< GeneralGA > &GA, const MTAnisoGAConf Configuration)
Definition: mtanisoga.cpp:42
This is a special objective function to fit phase tensor MT data.
Implements binary tournament selection for genetic algorithms.
void f(vector< double > &v1, vector< double > &v2, vector< double > &v3, vector< double > &v4)
Definition: perftest.cpp:17
A population that is encoded as a simple binary string.
AnnealingGA implements a genetic algorithm with an annealing style objective function.
Definition: AnnealingGA.h:15
virtual void GetData(const std::string filename)
read in data from file, determines format by ending
Definition: MTStation.cpp:597
void h(vector< double > &v1, vector< double > &v2, vector< double > &v3, vector< double > &v4)
Definition: perftest.cpp:40
void SetParams(const double mutation, const double crossover)
boost::shared_ptr< PlottableObjective > MTObjective
Definition: cadianiso.cpp:16
virtual ttranscribed GetValues(const tpopmember &member)
Re-Implements the function from BinaryTranscribe.
std::vector< gplib::rvec > tparamvector
We provide some typedefs that are used in other parts as well.
Definition: GeneralGA.h:29
This is the standard propagation class that generates a new population from the old one...
The class MTStation is used to store the transfer functions and related information for a MT-site...
Definition: MTStation.h:17
ublas::vector< bool > tpopmember
Definition: gentypes.h:22
std::vector< std::vector< int > > tparamindv
Definition: GeneralGA.h:32
int main(int argc, char *argv[])
Program to invert MT data for 1D anisotropic structure with a genetic algorithm.
Definition: mtanisoga.cpp:62
CMTStation MTData
Definition: cadijoint.cpp:15
void GetData(const std::string &filename)
void PrintAll(std::ostream &output)
Definition: UniquePop.cpp:36
This class stores a single unique copy of each population member that is added to it...
Definition: UniquePop.h:17
tgatype
Definition: mtanisoga.cpp:34
This class implements the Gray code representation of a binary string and the corresponding transcrip...
string version
Definition: makeinput.cpp:16
ublas::vector< double > tfitvec
Definition: gentypes.h:23
Implements a genetic algorithm based on the concept of pareto-optimality, best suited for multi-objec...
Definition: ParetoGA.h:17
CLevanisoConf Configuration
Definition: cadianiso.cpp:17
ublas::vector< int > tsizev
Definition: gentypes.h:27
void SetParams(const double InitT, const int AnnealG, const double AnnealR)
Set the parameters for the annealing process.
Definition: AnnealingGA.cpp:35
The basic exception class for all errors that arise in gplib.