ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/ZConsWriter.cpp
Revision: 1930
Committed: Wed Jan 12 22:41:40 2005 UTC (19 years, 5 months ago) by gezelter
File size: 6099 byte(s)
Log Message:
merging new_design branch into OOPSE-2.0

File Contents

# Content
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
50
51 namespace oopse {
52 ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) : info_(info) {
53 //use master - slave mode, only master node writes to disk
54 #ifdef IS_MPI
55 if(worldRank == 0){
56 #endif
57
58 output_.open(filename.c_str());
59
60 if(!output_){
61 sprintf( painCave.errMsg,
62 "Could not open %s for z constrain output_ \n", filename.c_str());
63 painCave.isFatal = 1;
64 simError();
65 }
66
67 output_ << "//time(fs)" << std::endl;
68 output_ << "//number of fixed z-constrain molecules" << std::endl;
69 output_ << "//global Index of molecule\tzconstrain force\tcurrentZPos" << std::endl;
70
71 #ifdef IS_MPI
72 }
73 #endif
74
75 }
76
77 ZConsWriter::~ZConsWriter()
78 {
79
80 #ifdef IS_MPI
81 if(worldRank == 0 ){
82 #endif
83 output_.close();
84 #ifdef IS_MPI
85 }
86 #endif
87 }
88
89 void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols){
90 #ifndef IS_MPI
91 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl;
92 output_ << fixedZmols.size() << std::endl;
93
94 std::list<ZconstraintMol>::const_iterator i;
95 for ( i = fixedZmols.begin(); i != fixedZmols.end(); ++i) {
96 output_ << i->mol->getGlobalIndex() <<"\t" << i->fz << "\t" << i->zpos << "\t" << i->param.zTargetPos <<std::endl;
97 }
98 #else
99 int nproc;
100 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
101 const int masterNode = 0;
102 int myNode = worldRank;
103 std::vector<int> tmpNFixedZmols(nproc, 0);
104 std::vector<int> nFixedZmolsInProc(nproc, 0);
105 tmpNFixedZmols[myNode] = fixedZmols.size();
106
107 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
108 MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0], nproc, MPI_INT,
109 MPI_SUM, MPI_COMM_WORLD);
110
111 MPI_Status ierr;
112 int zmolIndex;
113 double data[3];
114
115 if (masterNode == 0) {
116
117 std::vector<ZconsData> zconsData;
118 ZconsData tmpData;
119 for(int i =0 ; i < nproc; ++i) {
120 if (i == masterNode) {
121 std::list<ZconstraintMol>::const_iterator j;
122 for ( j = fixedZmols.begin(); j != fixedZmols.end(); ++j) {
123 tmpData.zmolIndex = j->mol->getGlobalIndex() ;
124 tmpData.zforce= j->fz;
125 tmpData.zpos = j->zpos;
126 tmpData.zconsPos = j->param.zTargetPos;
127 zconsData.push_back(tmpData);
128 }
129
130 } else {
131 for(int k =0 ; k < nFixedZmolsInProc[i]; ++k) {
132 MPI_Recv(&zmolIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD,&ierr);
133 MPI_Recv(data, 3, MPI_DOUBLE, i, 0, MPI_COMM_WORLD,&ierr);
134 tmpData.zmolIndex = zmolIndex;
135 tmpData.zforce= data[0];
136 tmpData.zpos = data[1];
137 tmpData.zconsPos = data[2];
138 zconsData.push_back(tmpData);
139 }
140 }
141
142 }
143
144
145 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl;
146 output_ << zconsData.size() << std::endl;
147
148 std::vector<ZconsData>::iterator l;
149 for (l = zconsData.begin(); l != zconsData.end(); ++l) {
150 output_ << l->zmolIndex << "\t" << l->zforce << "\t" << l->zpos << "\t" << l->zconsPos << std::endl;
151 }
152
153 } else {
154
155 std::list<ZconstraintMol>::const_iterator j;
156 for (j = fixedZmols.begin(); j != fixedZmols.end(); ++j) {
157 zmolIndex = j->mol->getGlobalIndex();
158 data[0] = j->fz;
159 data[1] = j->zpos;
160 data[2] = j->param.zTargetPos;
161 MPI_Send(&zmolIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
162 MPI_Send(data, 3, MPI_DOUBLE, masterNode, 0, MPI_COMM_WORLD);
163
164 }
165 }
166 #endif
167 }
168
169 }