ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/math/ParallelRandNumGen.cpp
Revision: 1793
Committed: Fri Aug 31 21:16:10 2012 UTC (12 years, 8 months ago) by gezelter
File size: 5063 byte(s)
Log Message:
Cleaning up some warning messages on linux

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     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (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 tim 392 MPI_Bcast(&seed, 1, MPI_UNSIGNED_LONG, masterNode, MPI_COMM_WORLD);
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 1313 MPI_Comm_size( MPI_COMM_WORLD, &nProcessors);
71 tim 398 MPI_Comm_rank( MPI_COMM_WORLD, &myRank_);
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     std::vector<uint32> bigSeed;
89     int nProcessors;
90 gezelter 1241 #ifdef IS_MPI
91 gezelter 1793 const int masterNode = 0;
92 gezelter 1313 MPI_Comm_size( MPI_COMM_WORLD, &nProcessors);
93 tim 398 MPI_Comm_rank( MPI_COMM_WORLD, &myRank_);
94 gezelter 1241 #else
95     nProcessors = 1;
96     myRank_ = 0;
97     #endif
98 tim 398 mtRand_ = new MTRand(nProcessors, myRank_);
99 tim 388
100 gezelter 1241 seed(); /** @todo calling virtual function in constructor is
101     not a good design */
102 gezelter 507 }
103 tim 388
104    
105 gezelter 507 void ParallelRandNumGen::seed( const uint32 oneSeed ) {
106 tim 388
107 gezelter 1313 unsigned long seed = oneSeed;
108 gezelter 1241 #ifdef IS_MPI
109 gezelter 1793 const int masterNode = 0;
110 tim 388 MPI_Bcast(&seed, 1, MPI_UNSIGNED_LONG, masterNode, MPI_COMM_WORLD);
111 gezelter 1241 #endif
112 tim 388 if (seed != oneSeed) {
113 gezelter 507 sprintf(painCave.errMsg,
114     "Using different seed to initialize ParallelRandNumGen.\n");
115     painCave.isFatal = 1;;
116     simError();
117 tim 388 }
118 gezelter 1241
119 gezelter 1313 unsigned long newSeed = oneSeed +nCreatedRNG_;
120 tim 388 mtRand_->seed(newSeed);
121 gezelter 1241
122 tim 388 ++nCreatedRNG_;
123 gezelter 507 }
124 tim 388
125 gezelter 507 void ParallelRandNumGen::seed() {
126 tim 388
127     std::vector<uint32> bigSeed;
128     int size;
129 gezelter 1793
130     #ifdef IS_MPI
131 tim 388 const int masterNode = 0;
132     if (worldRank == masterNode) {
133 gezelter 1241 #endif
134    
135 gezelter 507 bigSeed = mtRand_->generateSeeds();
136     size = bigSeed.size();
137 gezelter 1241
138     #ifdef IS_MPI
139 gezelter 507 MPI_Bcast(&size, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
140     MPI_Bcast(&bigSeed[0], size, MPI_UNSIGNED_LONG, masterNode, MPI_COMM_WORLD);
141 tim 388 }else {
142 gezelter 507 MPI_Bcast(&size, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
143     bigSeed.resize(size);
144     MPI_Bcast(&bigSeed[0], size, MPI_UNSIGNED_LONG, masterNode, MPI_COMM_WORLD);
145 tim 388 }
146 gezelter 1241 #endif
147 tim 388
148     if (bigSeed.size() == 1) {
149 gezelter 507 mtRand_->seed(bigSeed[0]);
150 tim 388 } else {
151 gezelter 507 mtRand_->seed(&bigSeed[0], bigSeed.size());
152 tim 388 }
153    
154     ++nCreatedRNG_;
155 gezelter 1241 }
156 tim 388 }

Properties

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