ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/ZConsWriter.cpp
(Generate patch)

Comparing trunk/src/io/ZConsWriter.cpp (file contents):
Revision 246 by gezelter, Wed Jan 12 22:41:40 2005 UTC vs.
Revision 1981 by gezelter, Mon Apr 14 18:32:51 2014 UTC

# Line 1 | Line 1
1 < /*
1 > /*
2   * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3   *
4   * The University of Notre Dame grants you ("Licensee") a
# Line 6 | Line 6
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
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
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.
# Line 37 | Line 28
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  
46
51   #include "io/ZConsWriter.hpp"
52   #include "utils/simError.h"
53  
54 <
55 < namespace oopse {
56 < ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) : info_(info) {
53 <  //use master - slave mode, only master node writes to disk
54 > namespace OpenMD {
55 >  ZConsWriter::ZConsWriter(SimInfo* info, const std::string& filename) : info_(info) {
56 >    //use master - slave mode, only master node writes to disk
57   #ifdef IS_MPI
58      if(worldRank == 0){
59   #endif
60  
61 <    output_.open(filename.c_str());
61 >      output_.open(filename.c_str());
62  
63 <    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 <    }
63 >      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  
70 <    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;
70 >      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  
74   #ifdef IS_MPI
75      }
76   #endif  
77  
78 < }
78 >  }
79  
80 < ZConsWriter::~ZConsWriter()
81 < {
80 >  ZConsWriter::~ZConsWriter()
81 >  {
82  
83   #ifdef IS_MPI
84 <  if(worldRank == 0 ){
84 >    if(worldRank == 0 ){
85   #endif  
86 <  output_.close();  
86 >      output_.close();  
87   #ifdef IS_MPI  
88 <  }
88 >    }
89   #endif
90 < }
90 >  }
91  
92 < void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols){
92 >  void ZConsWriter::writeFZ(const std::list<ZconstraintMol>& fixedZmols){
93   #ifndef IS_MPI
94      output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl;
95      output_ << fixedZmols.size() << std::endl;
96  
97      std::list<ZconstraintMol>::const_iterator i;
98      for ( i = fixedZmols.begin(); i != fixedZmols.end(); ++i) {
99 <        output_ << i->mol->getGlobalIndex() <<"\t" << i->fz << "\t" << i->zpos << "\t" << i->param.zTargetPos <<std::endl;
99 >      output_ << i->mol->getGlobalIndex() <<"\t" << i->fz << "\t" << i->zpos << "\t" << i->param.zTargetPos <<std::endl;
100      }
101   #else
102 <    int nproc;
100 <    MPI_Comm_size(MPI_COMM_WORLD, &nproc);
102 >  
103      const int masterNode = 0;
104 <    int myNode = worldRank;
104 >    int nproc;
105 >    int myNode;      
106 >    MPI_Comm_size( MPI_COMM_WORLD, &nproc);
107 >    MPI_Comm_rank( MPI_COMM_WORLD, &myNode);
108 >
109      std::vector<int> tmpNFixedZmols(nproc, 0);
110      std::vector<int> nFixedZmolsInProc(nproc, 0);
111      tmpNFixedZmols[myNode] = fixedZmols.size();
112      
113 <    //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
114 <    MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0], nproc, MPI_INT,
115 <                  MPI_SUM, MPI_COMM_WORLD);
113 >    //do MPI_ALLREDUCE to exchange the total number of atoms,
114 >    //rigidbodies and cutoff groups
115 >    MPI_Allreduce(&tmpNFixedZmols[0], &nFixedZmolsInProc[0],
116 >                  nproc, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
117  
118 <    MPI_Status ierr;
118 >    MPI_Status* ierr;
119      int zmolIndex;
120 <    double data[3];
120 >    RealType data[3];
121      
122 <    if (masterNode == 0) {
122 >    if (myNode == masterNode) {
123  
124 <        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 <                }                
124 >      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  
137 <            } else {
138 <                for(int k =0 ; k < nFixedZmolsInProc[i]; ++k) {
139 <                    MPI_Recv(&zmolIndex, 1, MPI_INT, i, 0, MPI_COMM_WORLD,&ierr);
140 <                    MPI_Recv(data, 3, MPI_DOUBLE, i, 0, MPI_COMM_WORLD,&ierr);
141 <                    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 <            }
137 >        } else {
138 >          for(int k =0 ; k < nFixedZmolsInProc[i]; ++k) {
139 >            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 >            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              
149 <        }
149 >      }
150  
151  
152 <        output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl;
153 <        output_ << zconsData.size() << std::endl;
152 >      output_ << info_->getSnapshotManager()->getCurrentSnapshot()->getTime() << std::endl;
153 >      output_ << zconsData.size() << std::endl;
154  
155 <        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 <        }
155 >      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          
160      } else {
161  
162 <        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 <            MPI_Send(&zmolIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
169 <            MPI_Send(data, 3, MPI_DOUBLE, masterNode, 0, MPI_COMM_WORLD);
162 >      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 >        MPI_Send(&zmolIndex, 1, MPI_INT, masterNode, 0, MPI_COMM_WORLD);
169 >        MPI_Send(data, 3, MPI_REALTYPE, masterNode, 0, MPI_COMM_WORLD);
170              
171 <        }
171 >      }
172      }
173   #endif
174 < }
174 >  }
175  
176   }

Comparing trunk/src/io/ZConsWriter.cpp (property svn:keywords):
Revision 246 by gezelter, Wed Jan 12 22:41:40 2005 UTC vs.
Revision 1981 by gezelter, Mon Apr 14 18:32:51 2014 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines