ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/ZConsWriter.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/ZConsWriter.cpp (file contents):
Revision 699 by tim, Fri Aug 15 19:24:13 2003 UTC vs.
Revision 1091 by tim, Tue Mar 16 19:22:56 2004 UTC

# Line 18 | Line 18 | ZConsWriter::ZConsWriter(const char* filename, vector<
18    
19     if(!output){
20       sprintf( painCave.errMsg,
21 <              "Could not open \"s\" for z constrain output \n",
22 <               filename);
21 >              "Could not open %s for z constrain output \n",
22 >         filename);
23       painCave.isFatal = 1;
24       simError();
25     }
26     output << "#number of z constrain molecules" << endl;
27     output << "#global Index of molecule\tzPos" << endl;
28 <
29 <   parameters = thePara;
30 <        writeZPos();
31 <
28 >   output << "#every frame will contain below data" <<endl;
29     output << "#time(fs)" << endl;
30     output << "#number of fixed z-constrain molecules" << endl;
31     output << "#global Index of molecule\tzconstrain force\tcurrentZPos" << endl;
32  
33 +   parameters = thePara;
34 +   writeZPos();
35  
36   #ifdef IS_MPI
37    }
# Line 52 | Line 51 | void ZConsWriter::writeFZ(double time, int num, int* i
51   #endif
52   }
53  
54 < void ZConsWriter::writeFZ(double time, int num, int* index, double* fz, double* curZPos)
55 < {
54 > /**
55 > *
56 > */
57 > void ZConsWriter::writeFZ(double time, int num, int* index, double* fz, double* curZPos, double* zpos){
58  
59   #ifndef IS_MPI
60    output << time << endl;
61    output << num << endl;
62 <        
62 >  
63    for(int i = 0; i < num; i++)
64 <    output << index[i] <<"\t" << fz[i] << "\t" << curZPos[i] << endl;
64 >    output << index[i] <<"\t" << fz[i] << "\t" << curZPos[i] << "\t" << zpos[i] <<endl;
65  
66   #else
67 <
67 >  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    int whichNode;
76    enum CommType { RequesPosAndForce, EndOfRequest} status;
77    double pos;
78    double force;
79 +  double zconsPos;
80    int localIndex;
81    MPI_Status ierr;
82    int tag = 0;
# Line 78 | Line 87 | void ZConsWriter::writeFZ(double time, int num, int* i
87      int *MolToProcMap;
88      MolToProcMap = mpiSim->getMolToProcMap();
89      
90 <    for(int i = 0; i < parameters->size(); i++){
90 >    for(int i = 0; i < (int)(parameters->size()); i++){
91        
92        globalIndexOfCurMol = (*parameters)[i].zconsIndex;
93        whichNode = MolToProcMap[globalIndexOfCurMol];
94        
95        if(whichNode == 0){
96          
97 <             for(int j = 0; j < num; j++)
98 <              if(index[j] == globalIndexOfCurMol){
99 <                localIndex = j;
100 <                break;
101 <              }
97 >       for(int j = 0; j < num; j++)
98 >        if(index[j] == globalIndexOfCurMol){
99 >          localIndex = j;
100 >          break;
101 >        }
102  
103 <                        force = fz[localIndex];
104 <                        pos = curZPos[localIndex];
105 <                  
103 >      force = fz[localIndex];
104 >      pos = curZPos[localIndex];
105 >      
106        }
107        else{
108          status = RequesPosAndForce;
109          MPI_Send(&status, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
110 <             MPI_Send(&globalIndexOfCurMol, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
111 <        MPI_Recv(&force, 1, MPI_DOUBLE_PRECISION, whichNode, tag, MPI_COMM_WORLD, &ierr);
112 <        MPI_Recv(&pos, 1, MPI_DOUBLE_PRECISION, whichNode, tag, MPI_COMM_WORLD, &ierr);
110 >        MPI_Send(&globalIndexOfCurMol, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD);
111 >        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 >        MPI_Recv(&zconsPos, 1, MPI_DOUBLE, whichNode, tag, MPI_COMM_WORLD, &ierr);
114        }
115  
116 <                output << globalIndexOfCurMol << "\t" << force << "\t" << pos << endl;
116 >     output << globalIndexOfCurMol << "\t" << force << "\t" << pos << "\t"<<  zconsPos << endl;
117                
118      } //End of Request Loop
119      
# Line 124 | Line 134 | void ZConsWriter::writeFZ(double time, int num, int* i
134      
135        switch (status){
136            
137 <        case RequesPosAndForce :
137 >         case RequesPosAndForce :
138            
139 <               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 <               MPI_Send(&fz[localIndex], 1, MPI_DOUBLE_PRECISION, 0, tag, MPI_COMM_WORLD);        
148 <               MPI_Send(&curZPos[localIndex], 1, MPI_DOUBLE_PRECISION, 0, tag, MPI_COMM_WORLD);          
149 <               break;
150 <            
139 >           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 >           MPI_Send(&fz[localIndex], 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);    
148 >           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 >           break;
151 >      
152          case EndOfRequest :
153          
154 <               done = true;
155 <               break;
154 >         done = true;
155 >         break;
156        }
157        
158      }
# Line 152 | Line 163 | void ZConsWriter::writeZPos(){
163  
164   }
165  
166 + /*
167 + *
168 + */
169   void ZConsWriter::writeZPos(){
170  
171   #ifdef IS_MPI
172    if(worldRank == 0){
173   #endif
174      
175 <     output << parameters->size() << endl;    
175 >    output << parameters->size() << endl;    
176      
177 <    for(int i =0 ; i < parameters->size(); i++)
177 >    for(int i =0 ; i < (int)(parameters->size()); i++)
178        output << (*parameters)[i].zconsIndex << "\t" <<  (*parameters)[i].zPos << endl;
179  
180   #ifdef IS_MPI

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines