libatomprobe
Library for Atom Probe Tomography (APT) computation
sampling.cpp
Go to the documentation of this file.
1 /*
2  * sampling.cpp - statistical sampling functions
3  * Copyright (C) 2015, D Haley
4 
5  * This program is free software: you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation, either version 3 of the License, or
8  * (at your option) any later version.
9 
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14 
15  * You should have received a copy of the GNU General Public License
16  * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18 
20 
22 
23 #include "gsl/gsl_randist.h"
24 
25 
26 namespace AtomProbe
27 {
28 
29 using std::vector;
30 void sampleIons(const vector<IonHit> &ions, float sampleFactor,
31  vector<IonHit> &sampled, bool strongRandom)
32 {
33  //obtain library random number generator, if we will use it
34  // otherwise, sampling function will use its own weak generator
35  gsl_rng *r;
36  if(strongRandom)
37  r=randGen.getRng();
38  else
39  r=0;
40 
41  size_t nSamples =sampleFactor*ions.size();
42  sampled.resize(nSamples);
43  randomSelect(sampled,ions,nSamples,r);
44 }
45 
46 #ifdef DEBUG
47 bool testSampling()
48 {
49  vector<IonHit> ions,sampled;
50  IonHit h(Point3D(1.0f,2.0f,3.0f),1.0f);
51  for(unsigned int ui=0;ui<10;ui++)
52  ions.push_back(h);
53 
54  sampleIons(ions,0.5,sampled);
55 
56  TEST(sampled.size() == 5,"sample count");
57  for(unsigned int ui=0;ui<sampled.size();ui++)
58  {
59  TEST(sampled[ui] == h,"sampled data integrity");
60  }
61 
62 
63  //Try again with strong randomsiation
64  sampled.clear();
65  sampleIons(ions,0.5,sampled,true);
66  TEST(sampled.size() == 5,"sample count");
67  for(unsigned int ui=0;ui<sampled.size();ui++)
68  {
69  TEST(sampled[ui] == h,"sampled data integrity");
70  }
71 
72  return true;
73 }
74 #endif
75 
76 }
#define TEST(f, g)
Definition: aptAssert.h:49
A 3D point data class storage.
Definition: point3D.h:39
void sampleIons(const std::vector< IonHit > &ions, float sampleFactor, std::vector< IonHit > &sampled, bool strongRandom=true)
Definition: sampling.cpp:30
size_t randomSelect(std::vector< T > &result, const std::vector< T > &source, size_t num, gsl_rng *rng)
Definition: sampling.h:40
This is a data holding class for POS file ions, from.
Definition: ionHit.h:36
gsl_rng * getRng() const
Obtain a GSL random number generator.
Definition: atomprobe.h:115
RandNumGen randGen
Definition: atomprobe.cpp:29