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. Acknowledgement of the program authors must be made in any |
10 |
* publication of scientific results based in part on use of the |
11 |
* program. An acceptable form of acknowledgement is citation of |
12 |
* the article in which the program was described (Matthew |
13 |
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
* Parallel Simulation Engine for Molecular Dynamics," |
16 |
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
* |
18 |
* 2. Redistributions of source code must retain the above copyright |
19 |
* notice, this list of conditions and the following disclaimer. |
20 |
* |
21 |
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
* notice, this list of conditions and the following disclaimer in the |
23 |
* documentation and/or other materials provided with the |
24 |
* distribution. |
25 |
* |
26 |
* This software is provided "AS IS," without a warranty of any |
27 |
* kind. All express or implied conditions, representations and |
28 |
* warranties, including any implied warranty of merchantability, |
29 |
* fitness for a particular purpose or non-infringement, are hereby |
30 |
* excluded. The University of Notre Dame and its licensors shall not |
31 |
* be liable for any damages suffered by licensee as a result of |
32 |
* using, modifying or distributing the software or its |
33 |
* derivatives. In no event will the University of Notre Dame or its |
34 |
* licensors be liable for any lost revenue, profit or data, or for |
35 |
* direct, indirect, special, consequential, incidental or punitive |
36 |
* damages, however caused and regardless of the theory of liability, |
37 |
* arising out of the use of or inability to use software, even if the |
38 |
* University of Notre Dame has been advised of the possibility of |
39 |
* such damages. |
40 |
*/ |
41 |
|
42 |
#include <algorithm> |
43 |
#include <iostream> |
44 |
#include <vector> |
45 |
|
46 |
|
47 |
#include "io/ZConsWriter.hpp" |
48 |
#include "utils/simError.h" |
49 |
#ifdef IS_MPI |
50 |
#include <mpi.h> |
51 |
#endif |
52 |
|
53 |
namespace oopse { |
54 |
ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) : info_(info) { |
55 |
//use master - slave mode, only master node writes to disk |
56 |
#ifdef IS_MPI |
57 |
if(worldRank == 0){ |
58 |
#endif |
59 |
|
60 |
output_.open(filename.c_str()); |
61 |
|
62 |
if(!output_){ |
63 |
sprintf( painCave.errMsg, |
64 |
"Could not open %s for z constrain output_ \n", filename.c_str()); |
65 |
painCave.isFatal = 1; |
66 |
simError(); |
67 |
} |
68 |
|
69 |
output_ << "//time(fs)" << std::endl; |
70 |
output_ << "//number of fixed z-constrain molecules" << std::endl; |
71 |
output_ << "//global Index of molecule\tzconstrain force\tcurrentZPos" << std::endl; |
72 |
|
73 |
#ifdef IS_MPI |
74 |
} |
75 |
#endif |
76 |
|
77 |
} |
78 |
|
79 |
ZConsWriter::~ZConsWriter() |
80 |
{ |
81 |
|
82 |
#ifdef IS_MPI |
83 |
if(worldRank == 0 ){ |
84 |
#endif |
85 |
output_.close(); |
86 |
#ifdef IS_MPI |
87 |
} |
88 |
#endif |
89 |
} |
90 |
|
91 |
void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols){ |
92 |
#ifndef IS_MPI |
93 |
output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl; |
94 |
output_ << fixedZmols.size() << std::endl; |
95 |
|
96 |
std::list<ZconstraintMol>::const_iterator i; |
97 |
for ( i = fixedZmols.begin(); i != fixedZmols.end(); ++i) { |
98 |
output_ << i->mol->getGlobalIndex() <<"\t" << i->fz << "\t" << i->zpos << "\t" << i->param.zTargetPos <<std::endl; |
99 |
} |
100 |
#else |
101 |
int nproc; |
102 |
MPI_Comm_size(MPI_COMM_WORLD, &nproc); |
103 |
const int masterNode = 0; |
104 |
int myNode = worldRank; |
105 |
std::vector<int> tmpNFixedZmols(nproc, 0); |
106 |
std::vector<int> nFixedZmolsInProc(nproc, 0); |
107 |
tmpNFixedZmols[myNode] = fixedZmols.size(); |
108 |
|
109 |
//do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups |
110 |
MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0], nproc, MPI_INT, |
111 |
MPI_SUM, MPI_COMM_WORLD); |
112 |
|
113 |
MPI_Status ierr; |
114 |
int zmolIndex; |
115 |
RealType data[3]; |
116 |
|
117 |
if (masterNode == 0) { |
118 |
|
119 |
std::vector<ZconsData> zconsData; |
120 |
ZconsData tmpData; |
121 |
for(int i =0 ; i < nproc; ++i) { |
122 |
if (i == masterNode) { |
123 |
std::list<ZconstraintMol>::const_iterator j; |
124 |
for ( j = fixedZmols.begin(); j != fixedZmols.end(); ++j) { |
125 |
tmpData.zmolIndex = j->mol->getGlobalIndex() ; |
126 |
tmpData.zforce= j->fz; |
127 |
tmpData.zpos = j->zpos; |
128 |
tmpData.zconsPos = j->param.zTargetPos; |
129 |
zconsData.push_back(tmpData); |
130 |
} |
131 |
|
132 |
} else { |
133 |
for(int k =0 ; k < nFixedZmolsInProc[i]; ++k) { |
134 |
MPI_Recv(&zmolIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD,&ierr); |
135 |
MPI_Recv(data, 3, MPI_REALTYPE, i, 0, MPI_COMM_WORLD,&ierr); |
136 |
tmpData.zmolIndex = zmolIndex; |
137 |
tmpData.zforce= data[0]; |
138 |
tmpData.zpos = data[1]; |
139 |
tmpData.zconsPos = data[2]; |
140 |
zconsData.push_back(tmpData); |
141 |
} |
142 |
} |
143 |
|
144 |
} |
145 |
|
146 |
|
147 |
output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl; |
148 |
output_ << zconsData.size() << std::endl; |
149 |
|
150 |
std::vector<ZconsData>::iterator l; |
151 |
for (l = zconsData.begin(); l != zconsData.end(); ++l) { |
152 |
output_ << l->zmolIndex << "\t" << l->zforce << "\t" << l->zpos << "\t" << l->zconsPos << std::endl; |
153 |
} |
154 |
|
155 |
} else { |
156 |
|
157 |
std::list<ZconstraintMol>::const_iterator j; |
158 |
for (j = fixedZmols.begin(); j != fixedZmols.end(); ++j) { |
159 |
zmolIndex = j->mol->getGlobalIndex(); |
160 |
data[0] = j->fz; |
161 |
data[1] = j->zpos; |
162 |
data[2] = j->param.zTargetPos; |
163 |
MPI_Send(&zmolIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD); |
164 |
MPI_Send(data, 3, MPI_REALTYPE, masterNode, 0, MPI_COMM_WORLD); |
165 |
|
166 |
} |
167 |
} |
168 |
#endif |
169 |
} |
170 |
|
171 |
} |