ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/ParallelRandNumGen.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 4949 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

File Contents

# User Rev Content
1 gezelter 507 /*
2 tim 388 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 tim 388 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 tim 388 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1782 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 tim 388 */
42    
43     #include "math/ParallelRandNumGen.hpp"
44     #ifdef IS_MPI
45     #include <mpi.h>
46 gezelter 1241 #endif
47 tim 388
48 gezelter 1390 namespace OpenMD {
49 tim 388
50 gezelter 507 int ParallelRandNumGen::nCreatedRNG_ = 0;
51 tim 388
52 gezelter 1241 ParallelRandNumGen::ParallelRandNumGen(const uint32& oneSeed) {
53 tim 388
54 gezelter 1313 unsigned long seed = oneSeed;
55    
56 gezelter 1241 #ifdef IS_MPI
57 gezelter 1793 const int masterNode = 0;
58 gezelter 1796 MPI::COMM_WORLD.Bcast(&seed, 1, MPI::UNSIGNED_LONG, masterNode);
59 gezelter 1241 #endif
60 gezelter 1313
61 tim 392 if (seed != oneSeed) {
62 gezelter 507 sprintf(painCave.errMsg,
63     "Using different seed to initialize ParallelRandNumGen.\n");
64     painCave.isFatal = 1;;
65     simError();
66 tim 392 }
67    
68 tim 388 int nProcessors;
69 gezelter 1241 #ifdef IS_MPI
70 gezelter 1796 nProcessors = MPI::COMM_WORLD.Get_size();
71     myRank_ = MPI::COMM_WORLD.Get_rank();
72 gezelter 1241 #else
73     nProcessors = 1;
74     myRank_ = 0;
75     #endif
76     //In order to generate independent random number stream, the
77     //actual seed used by random number generator is the seed passed
78     //to the constructor plus the number of random number generators
79     //which are already created.
80 gezelter 1313 unsigned long newSeed = oneSeed + nCreatedRNG_;
81 tim 398 mtRand_ = new MTRand(newSeed, nProcessors, myRank_);
82 gezelter 1241
83 tim 388 ++nCreatedRNG_;
84 gezelter 507 }
85 tim 388
86 gezelter 507 ParallelRandNumGen::ParallelRandNumGen() {
87 tim 388
88     int nProcessors;
89 gezelter 1241 #ifdef IS_MPI
90 gezelter 1796 nProcessors = MPI::COMM_WORLD.Get_size();
91     myRank_ = MPI::COMM_WORLD.Get_rank();
92 gezelter 1241 #else
93     nProcessors = 1;
94     myRank_ = 0;
95     #endif
96 tim 398 mtRand_ = new MTRand(nProcessors, myRank_);
97 tim 388
98 gezelter 1241 seed(); /** @todo calling virtual function in constructor is
99     not a good design */
100 gezelter 507 }
101 tim 388
102    
103 gezelter 507 void ParallelRandNumGen::seed( const uint32 oneSeed ) {
104 tim 388
105 gezelter 1313 unsigned long seed = oneSeed;
106 gezelter 1241 #ifdef IS_MPI
107 gezelter 1793 const int masterNode = 0;
108 gezelter 1796 MPI::COMM_WORLD.Bcast(&seed, 1, MPI::UNSIGNED_LONG, masterNode);
109 gezelter 1241 #endif
110 tim 388 if (seed != oneSeed) {
111 gezelter 507 sprintf(painCave.errMsg,
112     "Using different seed to initialize ParallelRandNumGen.\n");
113     painCave.isFatal = 1;;
114     simError();
115 tim 388 }
116 gezelter 1241
117 gezelter 1313 unsigned long newSeed = oneSeed +nCreatedRNG_;
118 tim 388 mtRand_->seed(newSeed);
119 gezelter 1241
120 tim 388 ++nCreatedRNG_;
121 gezelter 507 }
122 tim 388
123 gezelter 507 void ParallelRandNumGen::seed() {
124 tim 388
125     std::vector<uint32> bigSeed;
126 gezelter 1793
127     #ifdef IS_MPI
128 gezelter 1797 int size;
129 tim 388 const int masterNode = 0;
130     if (worldRank == masterNode) {
131 gezelter 1241 #endif
132    
133 gezelter 507 bigSeed = mtRand_->generateSeeds();
134 gezelter 1241
135     #ifdef IS_MPI
136 gezelter 1797 size = bigSeed.size();
137 gezelter 1796 MPI::COMM_WORLD.Bcast(&size, 1, MPI::INT, masterNode);
138     MPI::COMM_WORLD.Bcast(&bigSeed[0], size, MPI::UNSIGNED_LONG, masterNode);
139 tim 388 }else {
140 gezelter 1796 MPI::COMM_WORLD.Bcast(&size, 1, MPI::INT, masterNode);
141 gezelter 507 bigSeed.resize(size);
142 gezelter 1796 MPI::COMM_WORLD.Bcast(&bigSeed[0], size, MPI::UNSIGNED_LONG, masterNode);
143 tim 388 }
144 gezelter 1241 #endif
145 tim 388
146     if (bigSeed.size() == 1) {
147 gezelter 507 mtRand_->seed(bigSeed[0]);
148 tim 388 } else {
149 gezelter 507 mtRand_->seed(&bigSeed[0], bigSeed.size());
150 tim 388 }
151    
152     ++nCreatedRNG_;
153 gezelter 1241 }
154 tim 388 }

Properties

Name Value
svn:executable *
svn:keywords Author Id Revision Date