ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/ZConsWriter.cpp
Revision: 701
Committed: Wed Aug 20 14:34:04 2003 UTC (20 years, 10 months ago) by tim
File size: 4091 byte(s)
Log Message:
reformmating ZConstraint and fixe bug of error msg printing

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     "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 701 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     void ZConsWriter::writeFZ(double time, int num, int* index, double* fz, double* curZPos){
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     output << index[i] <<"\t" << fz[i] << "\t" << curZPos[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     int localIndex;
80     MPI_Status ierr;
81     int tag = 0;
82 tim 658
83 tim 699 if(worldRank == 0){
84 tim 658
85 tim 699 int globalIndexOfCurMol;
86     int *MolToProcMap;
87     MolToProcMap = mpiSim->getMolToProcMap();
88 tim 658
89 tim 699 for(int i = 0; i < parameters->size(); i++){
90 tim 658
91 tim 699 globalIndexOfCurMol = (*parameters)[i].zconsIndex;
92     whichNode = MolToProcMap[globalIndexOfCurMol];
93 tim 658
94 tim 699 if(whichNode == 0){
95 tim 658
96 tim 701 for(int j = 0; j < num; j++)
97     if(index[j] == globalIndexOfCurMol){
98     localIndex = j;
99     break;
100     }
101 tim 658
102 tim 701 force = fz[localIndex];
103     pos = curZPos[localIndex];
104    
105 tim 699 }
106     else{
107     status = RequesPosAndForce;
108     MPI_Send(&status, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
109 tim 701 MPI_Send(&globalIndexOfCurMol, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
110 tim 699 MPI_Recv(&force, 1, MPI_DOUBLE_PRECISION, whichNode, tag, MPI_COMM_WORLD, &ierr);
111     MPI_Recv(&pos, 1, MPI_DOUBLE_PRECISION, whichNode, tag, MPI_COMM_WORLD, &ierr);
112     }
113 tim 658
114 tim 701 output << globalIndexOfCurMol << "\t" << force << "\t" << pos << endl;
115 tim 699
116     } //End of Request Loop
117 tim 658
118 tim 699 //Send ending request message to slave nodes
119     status = EndOfRequest;
120     for(int i =1; i < mpiSim->getNumberProcessors(); i++)
121     MPI_Send(&status, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
122    
123     }
124     else{
125    
126     int whichMol;
127     bool done = false;
128    
129     while (!done){
130    
131     MPI_Recv(&status, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &ierr);
132 tim 658
133 tim 699 switch (status){
134    
135 tim 701 case RequesPosAndForce :
136 tim 699
137 tim 701 MPI_Recv(&whichMol, 1, MPI_INT, 0, tag, MPI_COMM_WORLD,&ierr);
138    
139     for(int i = 0; i < num; i++)
140     if(index[i] == whichMol){
141     localIndex = i;
142     break;
143     }
144    
145     MPI_Send(&fz[localIndex], 1, MPI_DOUBLE_PRECISION, 0, tag, MPI_COMM_WORLD);
146     MPI_Send(&curZPos[localIndex], 1, MPI_DOUBLE_PRECISION, 0, tag, MPI_COMM_WORLD);
147     break;
148    
149 tim 699 case EndOfRequest :
150    
151 tim 701 done = true;
152     break;
153 tim 699 }
154    
155 tim 658 }
156 tim 699
157 tim 658 }
158    
159     #endif
160    
161     }
162    
163 tim 701 /*
164     *
165     */
166 tim 699 void ZConsWriter::writeZPos(){
167 tim 658
168     #ifdef IS_MPI
169     if(worldRank == 0){
170     #endif
171    
172 tim 701 output << parameters->size() << endl;
173 tim 658
174 tim 699 for(int i =0 ; i < parameters->size(); i++)
175     output << (*parameters)[i].zconsIndex << "\t" << (*parameters)[i].zPos << endl;
176 tim 658
177     #ifdef IS_MPI
178     }
179     #endif
180     }