ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/ZConsWriter.cpp
Revision: 1198
Committed: Thu May 27 00:48:12 2004 UTC (20 years, 1 month ago) by tim
File size: 4293 byte(s)
Log Message:
in the progress of fixing MPI version of cutoff group

File Contents

# User Rev Content
1 tim 658 #include <algorithm>
2     #include <iostream>
3     #include <vector>
4     //#include <pair>
5     #include "ZConsWriter.hpp"
6     #include "simError.h"
7    
8     using namespace std;
9    
10 tim 699 ZConsWriter::ZConsWriter(const char* filename, vector<ZConsParaItem>* thePara)
11 tim 658 {
12     //use master - slave mode, only master node writes to disk
13     #ifdef IS_MPI
14     if(worldRank == 0){
15     #endif
16    
17     output.open(filename);
18    
19     if(!output){
20     sprintf( painCave.errMsg,
21 tim 738 "Could not open %s for z constrain output \n",
22 tim 701 filename);
23 tim 658 painCave.isFatal = 1;
24     simError();
25     }
26     output << "#number of z constrain molecules" << endl;
27 tim 699 output << "#global Index of molecule\tzPos" << endl;
28 tim 702 output << "#every frame will contain below data" <<endl;
29 tim 658 output << "#time(fs)" << endl;
30 tim 699 output << "#number of fixed z-constrain molecules" << endl;
31     output << "#global Index of molecule\tzconstrain force\tcurrentZPos" << endl;
32    
33 tim 701 parameters = thePara;
34     writeZPos();
35 tim 699
36 tim 658 #ifdef IS_MPI
37     }
38     #endif
39    
40     }
41    
42     ZConsWriter::~ZConsWriter()
43     {
44    
45     #ifdef IS_MPI
46     if(worldRank == 0 ){
47     #endif
48     output.close();
49     #ifdef IS_MPI
50     }
51     #endif
52     }
53    
54 tim 701 /**
55     *
56     */
57 tim 1091 void ZConsWriter::writeFZ(double time, int num, int* index, double* fz, double* curZPos, double* zpos){
58 tim 658
59     #ifndef IS_MPI
60     output << time << endl;
61     output << num << endl;
62 tim 701
63 tim 699 for(int i = 0; i < num; i++)
64 tim 1091 output << index[i] <<"\t" << fz[i] << "\t" << curZPos[i] << "\t" << zpos[i] <<endl;
65 tim 658
66 tim 699 #else
67 tim 701 int totalNum;
68     MPI_Allreduce(&num, &totalNum, 1, MPI_INT,MPI_SUM, MPI_COMM_WORLD);
69    
70     if(worldRank == 0){
71     output << time << endl;
72     output << totalNum << endl;
73     }
74    
75 tim 699 int whichNode;
76     enum CommType { RequesPosAndForce, EndOfRequest} status;
77     double pos;
78     double force;
79 tim 1091 double zconsPos;
80 tim 699 int localIndex;
81     MPI_Status ierr;
82     int tag = 0;
83 tim 658
84 tim 699 if(worldRank == 0){
85 tim 658
86 tim 699 int globalIndexOfCurMol;
87     int *MolToProcMap;
88     MolToProcMap = mpiSim->getMolToProcMap();
89 tim 658
90 mmeineke 787 for(int i = 0; i < (int)(parameters->size()); i++){
91 tim 658
92 tim 699 globalIndexOfCurMol = (*parameters)[i].zconsIndex;
93     whichNode = MolToProcMap[globalIndexOfCurMol];
94 tim 658
95 tim 699 if(whichNode == 0){
96 tim 658
97 tim 701 for(int j = 0; j < num; j++)
98     if(index[j] == globalIndexOfCurMol){
99     localIndex = j;
100     break;
101     }
102 tim 658
103 tim 701 force = fz[localIndex];
104     pos = curZPos[localIndex];
105    
106 tim 699 }
107     else{
108     status = RequesPosAndForce;
109     MPI_Send(&status, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
110 tim 701 MPI_Send(&globalIndexOfCurMol, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
111 gezelter 830 MPI_Recv(&force, 1, MPI_DOUBLE, whichNode, tag, MPI_COMM_WORLD, &ierr);
112     MPI_Recv(&pos, 1, MPI_DOUBLE, whichNode, tag, MPI_COMM_WORLD, &ierr);
113 tim 1091 MPI_Recv(&zconsPos, 1, MPI_DOUBLE, whichNode, tag, MPI_COMM_WORLD, &ierr);
114 tim 699 }
115 tim 658
116 tim 1091 output << globalIndexOfCurMol << "\t" << force << "\t" << pos << "\t"<< zconsPos << endl;
117 tim 699
118     } //End of Request Loop
119 tim 658
120 tim 699 //Send ending request message to slave nodes
121     status = EndOfRequest;
122 tim 1198 for(int i =1; i < mpiSim->getNprocessors(); i++)
123 tim 699 MPI_Send(&status, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
124    
125     }
126     else{
127    
128     int whichMol;
129     bool done = false;
130    
131     while (!done){
132    
133     MPI_Recv(&status, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &ierr);
134 tim 658
135 tim 699 switch (status){
136    
137 tim 701 case RequesPosAndForce :
138 tim 699
139 tim 701 MPI_Recv(&whichMol, 1, MPI_INT, 0, tag, MPI_COMM_WORLD,&ierr);
140    
141     for(int i = 0; i < num; i++)
142     if(index[i] == whichMol){
143     localIndex = i;
144     break;
145     }
146    
147 gezelter 830 MPI_Send(&fz[localIndex], 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
148 tim 1091 MPI_Send(&curZPos[localIndex], 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
149     MPI_Send(&zpos[localIndex], 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
150 tim 701 break;
151    
152 tim 699 case EndOfRequest :
153    
154 tim 701 done = true;
155     break;
156 tim 699 }
157    
158 tim 658 }
159 tim 699
160 tim 658 }
161    
162     #endif
163    
164     }
165    
166 tim 701 /*
167     *
168     */
169 tim 699 void ZConsWriter::writeZPos(){
170 tim 658
171     #ifdef IS_MPI
172     if(worldRank == 0){
173     #endif
174    
175 tim 701 output << parameters->size() << endl;
176 tim 658
177 mmeineke 787 for(int i =0 ; i < (int)(parameters->size()); i++)
178 tim 699 output << (*parameters)[i].zconsIndex << "\t" << (*parameters)[i].zPos << endl;
179 tim 658
180     #ifdef IS_MPI
181     }
182     #endif
183     }