ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/ConstraintWriter.cpp
Revision: 1983
Committed: Tue Apr 15 20:36:19 2014 UTC (11 years ago) by gezelter
File size: 6723 byte(s)
Log Message:
Fixes for ConstraintWriter in parallel, some cppcheck cleanup
starting preparation for 2.2 release

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. 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/ConstraintWriter.hpp"
52 #include "utils/simError.h"
53
54 // debugging:
55 #include <unistd.h>
56
57 namespace OpenMD {
58 ConstraintWriter::ConstraintWriter(SimInfo* info,
59 const std::string& filename): info_(info) {
60 //use master - slave mode, only master node writes to disk
61 #ifdef IS_MPI
62 if(worldRank == 0){
63 #endif
64 output_.open(filename.c_str());
65
66 if(!output_){
67 sprintf( painCave.errMsg,
68 "Could not open %s for Constraint output\n",
69 filename.c_str());
70 painCave.isFatal = 1;
71 simError();
72 }
73
74 output_ << "#time(fs)\t"
75 << "Index of atom 1\t"
76 << "Index of atom 2\tconstraint force" << std::endl;
77
78 #ifdef IS_MPI
79 }
80 #endif
81 }
82
83 ConstraintWriter::~ConstraintWriter() {
84 #ifdef IS_MPI
85 if(worldRank == 0 ){
86 #endif
87 output_.close();
88 #ifdef IS_MPI
89 }
90 #endif
91 }
92
93 void ConstraintWriter::writeConstraintForces(const std::list<ConstraintPair*>& constraints){
94 #ifndef IS_MPI
95 std::list<ConstraintPair*>::const_iterator i;
96 for ( i = constraints.begin(); i != constraints.end(); ++i) {
97
98 if ((*i)->getPrintForce()) {
99 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << "\t"
100 << (*i)->getConsElem1()->getGlobalIndex() <<"\t"
101 << (*i)->getConsElem2()->getGlobalIndex() <<"\t"
102 << (*i)->getConstraintForce() <<std::endl;
103 }
104 }
105 #else
106
107 const int masterNode = 0;
108 int nproc;
109 int myNode;
110 MPI_Comm_size( MPI_COMM_WORLD, &nproc);
111 MPI_Comm_rank( MPI_COMM_WORLD, &myNode);
112
113 std::vector<int> nConstraints(nproc, 0);
114 nConstraints[myNode] = constraints.size();
115
116 //do MPI_ALLREDUCE to exchange the total number of constraints:
117 MPI_Allreduce(MPI_IN_PLACE, &nConstraints[0], nproc, MPI_INT, MPI_SUM,
118 MPI_COMM_WORLD);
119
120 MPI_Status ierr;
121 int atom1, atom2, doPrint;
122 RealType force;
123
124 if (myNode == masterNode) {
125 std::vector<ConstraintData> constraintData;
126 ConstraintData tmpData;
127 for(int i = 0 ; i < nproc; ++i) {
128 if (i == masterNode) {
129 std::list<ConstraintPair*>::const_iterator j;
130 for ( j = constraints.begin(); j != constraints.end(); ++j) {
131 tmpData.atom1 = (*j)->getConsElem1()->getGlobalIndex();
132 tmpData.atom2 = (*j)->getConsElem2()->getGlobalIndex();
133 tmpData.constraintForce= (*j)->getConstraintForce();
134 tmpData.printForce = (*j)->getPrintForce();
135 constraintData.push_back(tmpData);
136 }
137
138 } else {
139 for(int k = 0; k < nConstraints[i]; ++k) {
140 MPI_Recv(&atom1, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
141 MPI_Recv(&atom2, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
142 MPI_Recv(&force, 1, MPI_REALTYPE, i, 0, MPI_COMM_WORLD, &ierr);
143 MPI_Recv(&doPrint, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &ierr);
144
145 tmpData.atom1 = atom1;
146 tmpData.atom2 = atom2;
147 tmpData.constraintForce = force;
148 tmpData.printForce = (doPrint == 0) ? false : true;
149 constraintData.push_back(tmpData);
150 }
151 }
152 }
153
154 std::vector<ConstraintData>::iterator l;
155 for (l = constraintData.begin(); l != constraintData.end(); ++l) {
156 if (l->printForce) {
157 output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << "\t"
158 << l->atom1 <<"\t"
159 << l->atom2 <<"\t"
160 << l->constraintForce << std::endl;
161 }
162 }
163 } else {
164
165 std::list<ConstraintPair*>::const_iterator j;
166 for (j = constraints.begin(); j != constraints.end(); ++j) {
167 int atom1 = (*j)->getConsElem1()->getGlobalIndex();
168 int atom2 = (*j)->getConsElem2()->getGlobalIndex();
169 RealType constraintForce= (*j)->getConstraintForce();
170 int printForce = (int) (*j)->getPrintForce();
171
172 MPI_Send(&atom1, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
173 MPI_Send(&atom2, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
174 MPI_Send(&constraintForce, 1, MPI_REALTYPE, masterNode, 0,
175 MPI_COMM_WORLD);
176 MPI_Send(&printForce, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
177 }
178 }
179 #endif
180 }
181 }