ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/ZConsWriter.cpp
Revision: 1969
Committed: Wed Feb 26 14:14:50 2014 UTC (11 years, 2 months ago) by gezelter
File size: 5947 byte(s)
Log Message:
Fixes to deal with deprecation of MPI C++ bindings.  We've reverted back to the
C calls.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date