1 |
/* |
2 |
* 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 |
* 1. Redistributions of source code must retain the above copyright |
10 |
* notice, this list of conditions and the following disclaimer. |
11 |
* |
12 |
* 2. Redistributions in binary form must reproduce the above copyright |
13 |
* 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 |
* |
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, 234107 (2008). |
39 |
* [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010). |
40 |
* [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). |
41 |
*/ |
42 |
|
43 |
#ifdef IS_MPI |
44 |
#include <mpi.h> |
45 |
#endif |
46 |
|
47 |
#include <algorithm> |
48 |
#include <iostream> |
49 |
#include <vector> |
50 |
|
51 |
#include "io/ZConsWriter.hpp" |
52 |
#include "utils/simError.h" |
53 |
|
54 |
namespace OpenMD { |
55 |
ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) : |
56 |
info_(info) { |
57 |
// use master - slave mode, only master node writes to disk |
58 |
#ifdef IS_MPI |
59 |
if(worldRank == 0){ |
60 |
#endif |
61 |
output_.open(filename.c_str()); |
62 |
|
63 |
if(!output_){ |
64 |
sprintf( painCave.errMsg, |
65 |
"Could not open %s for z constrain output_ \n", |
66 |
filename.c_str()); |
67 |
painCave.isFatal = 1; |
68 |
simError(); |
69 |
} |
70 |
|
71 |
output_ << "//time(fs)" << std::endl; |
72 |
output_ << "//number of fixed z-constrain molecules" << std::endl; |
73 |
output_ << "//global Index of molecule\tzconstrain force\tcurrentZPos" |
74 |
<< std::endl; |
75 |
#ifdef IS_MPI |
76 |
} |
77 |
#endif |
78 |
} |
79 |
|
80 |
ZConsWriter::~ZConsWriter() |
81 |
{ |
82 |
|
83 |
#ifdef IS_MPI |
84 |
if(worldRank == 0 ){ |
85 |
#endif |
86 |
output_.close(); |
87 |
#ifdef IS_MPI |
88 |
} |
89 |
#endif |
90 |
} |
91 |
|
92 |
void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols){ |
93 |
#ifndef IS_MPI |
94 |
output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() |
95 |
<< std::endl; |
96 |
output_ << fixedZmols.size() << std::endl; |
97 |
|
98 |
std::list<ZconstraintMol>::const_iterator i; |
99 |
for ( i = fixedZmols.begin(); i != fixedZmols.end(); ++i) { |
100 |
output_ << i->mol->getGlobalIndex() <<"\t" << i->fz << "\t" << i->zpos |
101 |
<< "\t" << i->param.zTargetPos <<std::endl; |
102 |
} |
103 |
#else |
104 |
|
105 |
const int masterNode = 0; |
106 |
int nproc; |
107 |
int myNode; |
108 |
MPI_Comm_size( MPI_COMM_WORLD, &nproc); |
109 |
MPI_Comm_rank( MPI_COMM_WORLD, &myNode); |
110 |
|
111 |
std::vector<int> tmpNFixedZmols(nproc, 0); |
112 |
std::vector<int> nFixedZmolsInProc(nproc, 0); |
113 |
tmpNFixedZmols[myNode] = fixedZmols.size(); |
114 |
|
115 |
//do MPI_ALLREDUCE to exchange the total number of atoms, |
116 |
//rigidbodies and cutoff groups |
117 |
MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0], |
118 |
nproc, MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
119 |
|
120 |
MPI_Status* ierr = NULL; |
121 |
int zmolIndex; |
122 |
RealType data[3]; |
123 |
|
124 |
if (myNode == masterNode) { |
125 |
|
126 |
std::vector<ZconsData> zconsData; |
127 |
ZconsData tmpData; |
128 |
for(int i =0 ; i < nproc; ++i) { |
129 |
if (i == masterNode) { |
130 |
std::list<ZconstraintMol>::const_iterator j; |
131 |
for ( j = fixedZmols.begin(); j != fixedZmols.end(); ++j) { |
132 |
tmpData.zmolIndex = j->mol->getGlobalIndex() ; |
133 |
tmpData.zforce= j->fz; |
134 |
tmpData.zpos = j->zpos; |
135 |
tmpData.zconsPos = j->param.zTargetPos; |
136 |
zconsData.push_back(tmpData); |
137 |
} |
138 |
|
139 |
} else { |
140 |
for(int k =0 ; k < nFixedZmolsInProc[i]; ++k) { |
141 |
MPI_Recv(&zmolIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD, ierr); |
142 |
MPI_Recv(data, 3, MPI_REALTYPE, i, 0, MPI_COMM_WORLD, ierr); |
143 |
tmpData.zmolIndex = zmolIndex; |
144 |
tmpData.zforce= data[0]; |
145 |
tmpData.zpos = data[1]; |
146 |
tmpData.zconsPos = data[2]; |
147 |
zconsData.push_back(tmpData); |
148 |
} |
149 |
} |
150 |
} |
151 |
|
152 |
output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() |
153 |
<< std::endl; |
154 |
output_ << zconsData.size() << std::endl; |
155 |
|
156 |
std::vector<ZconsData>::iterator l; |
157 |
for (l = zconsData.begin(); l != zconsData.end(); ++l) { |
158 |
output_ << l->zmolIndex << "\t" << l->zforce << "\t" << l->zpos |
159 |
<< "\t" << l->zconsPos << std::endl; |
160 |
} |
161 |
|
162 |
} else { |
163 |
|
164 |
std::list<ZconstraintMol>::const_iterator j; |
165 |
for (j = fixedZmols.begin(); j != fixedZmols.end(); ++j) { |
166 |
zmolIndex = j->mol->getGlobalIndex(); |
167 |
data[0] = j->fz; |
168 |
data[1] = j->zpos; |
169 |
data[2] = j->param.zTargetPos; |
170 |
MPI_Send(&zmolIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD); |
171 |
MPI_Send(data, 3, MPI_REALTYPE, masterNode, 0, MPI_COMM_WORLD); |
172 |
} |
173 |
} |
174 |
#endif |
175 |
} |
176 |
} |