ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/ConstraintWriter.cpp
Revision: 1984
Committed: Tue Apr 15 21:15:35 2014 UTC (11 years ago) by gezelter
File size: 6686 byte(s)
Log Message:
Some small compilation fixes for MSVC

File Contents

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