Main Page | Namespace List | Class Hierarchy | Class List | File List | Class Members | File Members

src/expde/evpar/communi.cc

Go to the documentation of this file.
00001 //    expde: expression templates for partial differential equations.
00002 //    Copyright (C) 2001  Christoph Pflaum
00003 //    This program is free software; you can redistribute it and/or modify
00004 //    it under the terms of the GNU General Public License as published by
00005 //    the Free Software Foundation; either version 2 of the License, or
00006 //    (at your option) any later version.
00007 //
00008 //    This program is distributed in the hope that it will be useful,
00009 //    but WITHOUT ANY WARRANTY; without even the implied warranty of
00010 //    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00011 //    GNU General Public License for more details.
00012 //
00013 //                 SEE  Notice1.doc made by 
00014 //                 LAWRENCE LIVERMORE NATIONAL LABORATORY
00015 //
00016 
00017 // ------------------------------------------------------------
00018 //
00019 // communi.cc
00020 //
00021 // ------------------------------------------------------------
00022 
00023 
00024 // Include level 0:
00025 
00026 #ifdef COMP_GNUOLD
00027  #include <iostream.h>
00028  #include <fstream.h>
00029  #include <time.h>
00030  #include <math.h>
00031 #else
00032  #include <iostream>
00033  #include <fstream>
00034  #include <ctime>
00035  #include <cmath>
00036 #endif 
00037 
00038 
00039 
00040 // Include level 0:
00041 #include "../parser.h"
00042 
00043 // Include level 1:
00044 #include "../paramete.h"
00045 #include "../abbrevi.h"
00046 #include "../math_lib/math_lib.h"
00047 
00048 // Include level 2:
00049 #include "../basic/basic.h"
00050 
00051 // Include level 3:
00052 #include "../domain/domain.h"
00053 
00054 // Include level 4:
00055 #include "../formulas/boundy.h"
00056 #include "../formulas/loc_sten.h"
00057 
00058 // Include level 5
00059 #include "../grid/gpar.h"
00060 #include "../grid/parallel.h"
00061 #include "../grid/mgcoeff.h"
00062 #include "../grid/sto_man.h"
00063 #include "../grid/gridbase.h"
00064 #include "../grid/grid.h"
00065 
00067 // communication of points
00069 
00070 
00071 /* old
00072 bool Grid_base::Send_I_point_in_direction(int d, 
00073                                           Point_hashtable1* poi, int t,
00074                                           bool also_for_prolongation,
00075                                           Index3D my_lev_index,
00076                                           Index3D next_index) {
00077   Index3D I;
00078   Index3D nextS;
00079 
00080   // new
00081   if(poi->typ >= interior   ||
00082      (poi->typ == multigrid && t<=poi->finest_level)) {
00083     I = poi->ind;
00084     if(t>=poi->Give_Tiefe()) {
00085       if(I << next_index) {
00086         return true;
00087       }
00088 
00089       if(d<6) {
00090         nextS = I.next((dir_3D)d,t);
00091       }
00092       else if(d<18) {
00093         nextS = I.next((Edges_cell)(d-6),t);
00094       }
00095       else
00096         nextS = I.next((dir_sons)(d-18),t);
00097       
00098       
00099       if(nextS << next_index) { // also on finest level for variable coefficients
00100         return true;
00101       }
00102       
00103       //      if(also_for_prolongation && t+1<=poi->my_finest_parallel_level) {
00104       if(also_for_prolongation) {
00105         if(d<6) {
00106           nextS = I.next((dir_3D)d,t+1);
00107         }
00108         else if(d<18) {
00109           nextS = I.next((Edges_cell)(d-6),t+1);
00110         }
00111         else
00112           nextS = I.next((dir_sons)(d-18),t+1);
00113         
00114         if(nextS < next_index) {
00115           return true;
00116         }
00117       }
00118     }
00119   }
00120   return false;
00121 }
00122 */
00123 
00124 
00125 bool Is_cell_in_direction(dir_sons c, dir_3D d) {
00126   if(d == Wdir || d == Edir) {
00127     if(xCoord(c) == d) 
00128       return true;
00129     else
00130       return false;
00131   }
00132   if(d == Ndir || d == Sdir) {
00133     if(yCoord(c) == d) 
00134       return true;
00135     else
00136       return false;
00137   }
00138   if(d == Tdir || d == Ddir) {
00139     if(zCoord(c) == d) 
00140       return true;
00141     else
00142       return false;
00143   }
00144   return false;
00145 }
00146 
00147 
00148 bool Grid_base::Send_B_point_in_direction(int d,
00149                                           Bo2Point* poi,
00150                                           Index3D next_index) {
00151   Index3D Icell;
00152   Index3D I, I_edge;
00153   int t;
00154   int ce, co;
00155 
00156   //  send one layer in left direction
00157   // this is for the following situation:
00158   //
00159   //     * 
00160   //     ^
00161   //  W  |  E
00162   //     +
00163   //
00164   // + is a boundary point from E which is needed in W
00165 
00166   I     = poi->ind;
00167   I_edge = I.next(poi->direction,max_level+1);
00168   t = max_level;
00169 
00170   // return false;
00171 
00172   // send two layers in right direction and one in left 
00173   if(I_edge < my_index) {
00174     // first part: for send on levels
00175     for(ce=0;ce<8;++ce) {
00176       if(Is_cell_in_direction((dir_sons)ce,poi->direction)) {
00177         // study relevant cell
00178         Icell = I.next((dir_sons)ce,t+1);
00179         for(co=0;co<8;++co) {
00180           // study corners of cell
00181           if(Icell.neighbour((dir_sons)co) << next_index)
00182             return true;
00183         }
00184       }
00185     }
00186     // second part: for restriction
00187     for(ce=0;ce<8;++ce) {
00188       if(Is_cell_in_direction((dir_sons)ce,poi->direction)) {
00189         // study relevant cell
00190         Icell = I.next((dir_sons)ce,t+1).father();
00191         for(co=0;co<8;++co) {
00192           // study corners of cell
00193           if(Icell.neighbour((dir_sons)co) << next_index)
00194             return true;
00195         }
00196       }
00197     }
00198   }
00199   return false;
00200 }
00201 
00202 
00203 bool Grid_base::Send_Z_point_in_direction(int d, BoCell* poi,
00204                                           Index3D next_index) {
00205   Index3D I, Icell;
00206   int co;
00207 
00208   if(poi->Exists_bocellpoint()) {
00209     I = poi->ind;
00210     if(I < my_index) {
00211       // first part: for send on levels
00212       Icell = I;
00213       for(co=0;co<8;++co) {
00214         // study corners of cell
00215         if(Icell.neighbour((dir_sons)co) << next_index)
00216           return true;
00217       }
00218       // second part: for restriction
00219       Icell = I.father();
00220       for(co=0;co<8;++co) {
00221         // study corners of cell
00222         if(Icell.neighbour((dir_sons)co) << next_index) {
00223           return true;
00224         }
00225       }
00226     }
00227   }
00228   return false;
00229 }
00230 
00231 
00232 void Grid_base::Prepare_communication() {
00233   // 1. initialize storage
00234   number_receive_face = new int*[max_level+1];
00235   number_send_face = new int*[max_level+1];
00236   number_receive_face_prolongation = new int*[max_level+1];
00237   number_send_face_prolongation = new int*[max_level+1];
00238 
00239   var_receive = new double***[max_level+1];
00240   var_send    = new double***[max_level+1];
00241   var_receive_prolongation = new double***[max_level+1];
00242   var_send_prolongation    = new double***[max_level+1];
00243 
00244   for(int i=0;i<max_level+1;++i) {
00245     number_receive_face[i] = new int[26];
00246     number_send_face[i] = new int[26];
00247     number_receive_face_prolongation[i] = new int[26];
00248     number_send_face_prolongation[i]    = new int[26];
00249 
00250     var_receive[i]              = new double**[26];
00251     var_send[i]                 = new double**[26];
00252     var_receive_prolongation[i] = new double**[26];
00253     var_send_prolongation[i]    = new double**[26];
00254 
00255     for(int k=0;k<26;++k) {
00256       number_receive_face[i][k] = 0;
00257       number_send_face[i][k] = 0;
00258       number_receive_face_prolongation[i][k] = 0;
00259       number_send_face_prolongation[i][k] = 0;
00260     }
00261   }
00262 
00263   // 2. prepare
00264   Prepare_communication_all_grids();
00265   Prepare_communication_coarser_grids();
00266 }
00267 
00268 
00269 bool Grid_base::Send_prol_point_in_direction(int d, 
00270                                              Point_hashtable1* poi, int t,
00271                                              Index3D my_lev_index,
00272                                              Index3D next_index) {
00273   Index3D I, nextS;
00274 
00275   // new
00276   // t point level of coarse grid points
00277   if( poi->typ >= interior   ||
00278      (poi->typ == multigrid && t<=poi->finest_level)) {
00279     I = poi->ind;
00280     if(t>=poi->Give_Tiefe()) {
00281       if(I << next_index) {
00282         return true;
00283       }
00284       if(d<6) {
00285         nextS = I.next((dir_3D)d,t+1);
00286       }
00287       else if(d<18) {
00288         nextS = I.next((Edges_cell)(d-6),t+1);
00289         }
00290       else
00291         nextS = I.next((dir_sons)(d-18),t+1);
00292       if(nextS << next_index) {
00293         return true;
00294       }
00295     }
00296   }
00297 
00298   return false;
00299 }
00300 
00301 // for prolongation operators on coarser grids
00302 void Grid_base::Prepare_communication_coarser_grids() {
00303   int d, num, l;
00304   Pointtype ptyp;
00305   MPI_Request req[2];
00306   MPI_Status status[2];
00307   int num_message;
00308   Index3D I;
00309   Index3D my_lev_index;
00310   Index3D next_index;
00311   
00312   int* send_info[26];
00313   int* receive_info[26];
00314 
00315   int num_rec, num_send;
00316 
00317   int rank_source, rank_destination;
00318 
00319   // l is level with respect to coarse points
00320   // send point level is 1 larger
00321   for(l=n_parallel-2;l>=min_level;--l) {
00322     for(d=0;d<26;++d) {
00323       // 0. Step:
00324       if(d<6) {
00325         rank_source      = give_next_rank_source((dir_3D)d,l+1);
00326         rank_destination = give_next_rank_destination((dir_3D)d,l+1);  
00327       }
00328       else if(d<18) {
00329         rank_source      = give_next_rank_source((Edges_cell)(d-6),l+1);
00330         rank_destination = give_next_rank_destination((Edges_cell)(d-6),l+1);  
00331       }
00332       else {
00333         rank_source      = give_next_rank_source((dir_sons)(d-18),l+1);
00334         rank_destination = give_next_rank_destination((dir_sons)(d-18),l+1);
00335       }
00336 
00337       // 1. Step: 
00338       my_lev_index = Give_my_level_index(l+2);
00339       //old
00340       //      next_index = Give_next_on_level(d,my_lev_index);
00341       //new
00342       next_index = Give_level_index_of(rank_destination,l+1);
00343       if(next_index == Index3D(3,3,3))
00344         rank_destination = -1;
00345       if(Give_my_coarsest_level()>=l+1)
00346         rank_source      = -1;
00347 
00348       // 2. Step: Calc number send for faces
00349       number_send_face_prolongation[l][d] = 0;
00350     
00351       // 2.1. Step: interior points
00352 
00353 
00354 
00355       if(l>=my_coarsest_level) {
00356         // allowed to send
00357         iterate_hash1 {
00358           if(Send_prol_point_in_direction(d,point1,l,
00359                                           my_lev_index,next_index)) {
00360             ++number_send_face_prolongation[l][d];
00361           }
00362         }
00363       }
00364 
00365 
00366 //number_send_face_prolongation[l][d] = 0;
00367 
00368 
00369 /*
00370   if(l==2 && d==22)
00371   cout << " my rank: " << my_rank
00372   << " rank_source: " << rank_source
00373   << " rank_destination: " << rank_destination
00374   << endl;
00375 */
00376 
00377       // 3. Step: Send number_faces
00378       num_rec  = 0;
00379       num_send = number_send_face_prolongation[l][d];
00380       num_message = 0;
00381       if(rank_source != -1) {
00382         MPI_Irecv(&num_rec,1,
00383                   MPI_INT,rank_source,200+d,comm,
00384                   &req[num_message]);
00385         ++num_message;
00386       }
00387       if(rank_destination != -1) {
00388         MPI_Isend(&num_send,1,
00389                   MPI_INT,rank_destination,200+d,comm,
00390                   &req[num_message]);
00391         ++num_message;
00392       }
00393       MPI_Waitall(num_message,req,status);
00394       number_receive_face_prolongation[l][d] = num_rec;
00395       
00396       // construct double array 
00397       if(number_send_face_prolongation[l][d]!=0) {
00398         var_send_prolongation[l][d] =
00399           new double*[number_send_face_prolongation[l][d]];
00400       }
00401       else {
00402         var_send_prolongation[l][d] = NULL;
00403       }
00404       if(number_receive_face_prolongation[l][d]!=0) {
00405         var_receive_prolongation[l][d] = 
00406           new double*[number_receive_face_prolongation[l][d]];
00407       }
00408       else {
00409         var_receive_prolongation[l][d] = NULL;
00410       }
00411       
00412     
00413       // 4. Step: info = new ...
00414       if(number_send_face_prolongation[l][d]>0) {
00415         send_info[d] = new int[3*number_send_face_prolongation[l][d]];
00416       }
00417       else {
00418         send_info[d] = NULL;
00419       }      
00420       if(number_receive_face_prolongation[l][d]>0) {
00421         receive_info[d] = new int[3*number_receive_face_prolongation[l][d]];
00422       }
00423       else {
00424         receive_info[d] = NULL;
00425       }
00426       
00427 
00428       
00429       // 5. Step: transmitter interior points
00430       // 5.1 build informations about points to be send
00431       // How informations are stored:
00432       // Index of point       
00433       // ind_x, ind_y, ind_z     >   3 int
00434       num = 0;
00435       if(l>=my_coarsest_level) {
00436         // allowed to send
00437         iterate_hash1 {
00438           if(Send_prol_point_in_direction(d,point1,l,
00439                                           my_lev_index,next_index)) {
00440             I = point1->ind;
00441             
00442             send_info[d][num*3]   = I.I_x().get();
00443             send_info[d][num*3+1] = I.I_y().get();
00444             send_info[d][num*3+2] = I.I_z().get();
00445 
00446             if(developer_version) {
00447               if(Give_variable_slow(I,l)==NULL) {
00448                 cout << " error in  5.1. Prepare_communication_coarser_grids()"
00449                      << endl;
00450               }
00451             }
00452 
00453             var_send_prolongation[l][d][num] = Give_variable(I,l);
00454             num++;
00455           }
00456         }
00457       }
00458 
00459         
00460       // 5.2 part: send points
00461       num_message = 0;
00462       if(rank_source != -1 && number_receive_face_prolongation[l][d]>0) {
00463         MPI_Irecv(receive_info[d] ,
00464                   3*number_receive_face_prolongation[l][d],
00465                   MPI_INT,rank_source,210+d,comm,
00466                   &req[num_message]);
00467         ++num_message;
00468       }
00469       if(rank_destination != -1 && number_send_face_prolongation[l][d]>0) {
00470         MPI_Isend(send_info[d],
00471                   3*number_send_face_prolongation[l][d],
00472                   MPI_INT,rank_destination,210+d,comm,
00473                   &req[num_message]);
00474         ++num_message;
00475       }
00476       MPI_Waitall(num_message,req,status);
00477  
00478          
00479       // 5.3 part: set points
00480       for(num=0;num<number_receive_face_prolongation[l][d];++num) {
00481         I = Index3D(receive_info[d][num*3],
00482                     receive_info[d][num*3+1],
00483                     receive_info[d][num*3+2]);
00484         var_receive_prolongation[l][d][num] = Give_variable_slow(I,l);
00485         if(var_receive_prolongation[l][d][num]==NULL) {
00486           if(Exists_Point(I)==false) {
00487             Add_point(I, parallel_p);
00488           }
00489           Add_double(I,l,Give_max_num_var());
00490           var_receive_prolongation[l][d][num] = Give_variable(I,l);
00491         }
00492 
00493 
00494 
00495         if(developer_version) {
00496           ptyp = Give_type(I);
00497           if(ptyp!=parallel_p && ptyp!=exterior) {
00498             cout << "\n error in 5.3 prolongation, my rank: " 
00499                  << my_rank 
00500                  << " Index: " 
00501                  << " x: " << receive_info[d][num*3]
00502                  << " y: " << receive_info[d][num*3+1]
00503                  << " z: " << receive_info[d][num*3+2]
00504                  << "typ: " << ptyp 
00505                  << endl;
00506           }
00507         }
00508       }
00509 
00510       // 6. part delete 
00511       if(send_info[d]!=NULL) {
00512         delete(send_info[d]);
00513         send_info[d] = NULL;
00514       }
00515       if(receive_info[d]!=NULL) {
00516         delete(receive_info[d]);
00517         receive_info[d] = NULL;
00518       }
00519     }
00520   }
00521 }
00522 
00523 
00524 
00525 
00526 bool Grid_base::Send_I_point_in_direction(int d, 
00527                                           Point_hashtable1* poi, int t,
00528                                           bool also_for_prolongation,
00529                                           Index3D my_lev_index,
00530                                           Index3D next_index) {
00531   Index3D I;
00532   Index3D nextS;
00533 
00534   I = poi->ind;
00535   if(I << next_index) {
00536     return true;
00537   }
00538   
00539   if(d<6) {
00540     nextS = I.next((dir_3D)d,t);
00541   }
00542   else if(d<18) {
00543     nextS = I.next((Edges_cell)(d-6),t);
00544   }
00545   else
00546     nextS = I.next((dir_sons)(d-18),t);
00547   
00548   
00549   if(nextS << next_index) { // also on finest level for variable coefficients
00550     return true;
00551   }
00552   
00553   //      if(also_for_prolongation && t+1<=poi->my_finest_parallel_level) {
00554   if(also_for_prolongation) {
00555     if(d<6) {
00556       nextS = I.next((dir_3D)d,t+1);
00557     }
00558     else if(d<18) {
00559       nextS = I.next((Edges_cell)(d-6),t+1);
00560     }
00561     else
00562       nextS = I.next((dir_sons)(d-18),t+1);
00563     
00564     if(nextS < next_index) {
00565       return true;
00566     }
00567   }
00568   return false;
00569 }
00570 
00571 
00572 
00573 void Grid_base::Prepare_communication_all_grids() {
00574   Index1D nextS;
00575   Pointtype ptyp;
00576   int d, num, l;
00577   MPI_Request req[12];
00578   MPI_Status status[12];
00579   int num_message;
00580   Index3D I, next_index;
00581   Index3D my_lev_index;
00582 
00583   int* send_info[26];
00584   int* receive_info[26];
00585   // I interior, B boundary, Z Zellfreiheitsgrad
00586   int num_rec_face_I[26],  num_rec_face_B[26],  num_rec_face_Z[26];
00587   int num_send_face_I[26], num_send_face_B[26], num_send_face_Z[26];
00588 
00589   int three_rec[26][3];
00590   int three_send[26][3];
00591 
00592   int rank_source, rank_destination;
00593   bool also_for_prolongation;
00594 
00595   for(l=Max_level();l>=my_coarsest_level;--l) {
00596   //for(l=Max_level();l>=Max_level();--l) {
00597 
00598     // -1. set my_lev_index 
00599     my_lev_index = Give_my_level_index(l+1);
00600     
00601 
00602     // construct extra array for invariant points which are needed.
00603     int lenght_invar_array = 0;
00604     iterate_hash1 {
00605       if(point1->typ >= interior   ||
00606          (point1->typ == multigrid && l<=point1->finest_level)) {
00607         if(l>=point1->Give_Tiefe()) {
00608           lenght_invar_array++;
00609         }
00610       }
00611     }
00612     Point_hashtable1** invariant_array = new Point_hashtable1*[lenght_invar_array];
00613     int array_i = 0;
00614     iterate_hash1 {
00615       if(point1->typ >= interior   ||
00616          (point1->typ == multigrid && l<=point1->finest_level)) {
00617         if(l>=point1->Give_Tiefe()) {
00618           invariant_array[array_i] = point1;
00619           array_i++;
00620         }
00621       }
00622     }
00623 
00624     for(d=0;d<26;++d) {
00625       // 0. Step: 
00626       next_index = Give_next_on_level(d,my_lev_index);
00627 
00628       // 1. Step:
00629       if(d<6) {
00630         rank_source      = give_next_rank_source((dir_3D)d,l);
00631         rank_destination = give_next_rank_destination((dir_3D)d,l);  
00632       }
00633       else if(d<18) {
00634         rank_source      = give_next_rank_source((Edges_cell)(d-6),l);
00635         rank_destination = give_next_rank_destination((Edges_cell)(d-6),l);  
00636       }
00637       else {
00638         rank_source      = give_next_rank_source((dir_sons)(d-18),l);
00639         rank_destination = give_next_rank_destination((dir_sons)(d-18),l);
00640       }
00641 
00642       // 2. Step: Calc number send for faces
00643       num_send_face_I[d] = 0;
00644       num_send_face_B[d] = 0;
00645       num_send_face_Z[d] = 0;
00646     
00647       // 2.1. Step: interior points
00648       // level_prolongation must be finer for 
00649       // prolongation operators on coarser grids
00650       if(l==Max_level() || l<n_parallel-1) {
00651         also_for_prolongation = false;
00652       }
00653       else {
00654         also_for_prolongation = true;
00655       }
00656 
00657       for(int i = 0;i < lenght_invar_array;i++) {
00658         if(Send_I_point_in_direction(d,invariant_array[i],l,
00659                                      also_for_prolongation,
00660                                      my_lev_index,
00661                                      next_index)) {
00662           ++num_send_face_I[d];
00663         }
00664       }
00665 
00666      /* old
00667       iterate_hash1 {
00668         if(Send_I_point_in_direction(d,point1,l,
00669                                      also_for_prolongation,
00670                                      my_lev_index,
00671                                      next_index)) {
00672           ++num_send_face_I[d];
00673         }
00674       }
00675      */
00676       if(l==Max_level()) {
00677         // 2.2. Step: boundary points
00678         iterate_hash3 {
00679           if(Send_B_point_in_direction((dir_3D)d,bo2point,next_index)) {
00680             ++num_send_face_B[d];
00681           }
00682         }
00683         // 2.3. Step: Zellfreiheitsgrad
00684         iterate_hash2 {
00685           if(Send_Z_point_in_direction(d,bocell,next_index)) {
00686             ++num_send_face_Z[d];
00687           }
00688         }
00689       }
00690       // 2.4. Step: sum up
00691       number_send_face[l][d] =
00692         num_send_face_I[d] +
00693         num_send_face_B[d] +
00694         num_send_face_Z[d];
00695       
00696       // 3. Step: Send number_faces
00697       num_message = 0;
00698       three_send[d][0]=num_send_face_I[d];
00699       three_send[d][1]=num_send_face_B[d];
00700       three_send[d][2]=num_send_face_Z[d];
00701       three_rec[d][0] =0;
00702       three_rec[d][1] =0;
00703       three_rec[d][2] =0;
00704       if(rank_source != -1) {
00705         MPI_Irecv(three_rec[d],3,
00706                   MPI_INT,rank_source,20+d,comm,
00707                   &req[num_message]);
00708         ++num_message;
00709       }
00710       if(rank_destination != -1) {
00711         MPI_Isend(three_send[d],3,
00712                   MPI_INT,rank_destination,20+d,comm,
00713                   &req[num_message]);
00714         ++num_message;
00715       }
00716       MPI_Waitall(num_message,req,status);
00717       
00718       num_rec_face_I[d]=three_rec[d][0];
00719       num_rec_face_B[d]=three_rec[d][1];
00720       num_rec_face_Z[d]=three_rec[d][2];
00721       
00722       // - sum up
00723       number_receive_face[l][d] =
00724         num_rec_face_I[d] +
00725         num_rec_face_B[d] +
00726         num_rec_face_Z[d];
00727       
00728       // construct double array 
00729       if(number_send_face[l][d]!=0)
00730         var_send[l][d] = new double*[number_send_face[l][d]];
00731       else var_send[l][d] = NULL;
00732       if(number_receive_face[l][d]!=0)
00733         var_receive[l][d] = 
00734           new double*[number_receive_face[l][d]];
00735       else var_receive[l][d] = NULL;
00736       
00737     
00738       // 4. Step: info = new ...
00739       if(MAX(3*num_send_face_I[d],
00740              4*num_send_face_B[d],
00741              3*num_send_face_Z[d])>0)
00742         send_info[d] = new int[MAX(3*num_send_face_I[d],
00743                                    4*num_send_face_B[d],
00744                                    3*num_send_face_Z[d])];
00745       else send_info[d]  = NULL;
00746       
00747       if(MAX(3*num_rec_face_I[d],
00748              4*num_rec_face_B[d],
00749              3*num_rec_face_Z[d])>0)
00750         receive_info[d] = new int[MAX(3*num_rec_face_I[d],
00751                                       4*num_rec_face_B[d],
00752                                       3*num_rec_face_Z[d])];
00753       else receive_info[d] = NULL;
00754       
00755       // 5. Step: transmitter interior points
00756       // 5.1 build informations about points to be send
00757       // How informations are stored:
00758       // Index of point       
00759       // ind_x, ind_y, ind_z     >   3 int
00760 
00761       num = 0;
00762       for(int i = 0;i < lenght_invar_array;i++) {
00763         point1 = invariant_array[i];
00764         if(Send_I_point_in_direction(d,point1,l,
00765                                      also_for_prolongation, my_lev_index,
00766                                      next_index)) {
00767           I = point1->ind;
00768           
00769           send_info[d][num*3]   = I.I_x().get();
00770           send_info[d][num*3+1] = I.I_y().get();
00771           send_info[d][num*3+2] = I.I_z().get();
00772           var_send[l][d][num] = Give_variable(I,l);
00773           num++;
00774         }
00775       }
00776 
00777       /* old
00778       num = 0;
00779       iterate_hash1 {
00780         if(Send_I_point_in_direction(d,point1,l,
00781                                      also_for_prolongation, my_lev_index,
00782                                      next_index)) {
00783           I = point1->ind;
00784           
00785           send_info[d][num*3]   = I.I_x().get();
00786           send_info[d][num*3+1] = I.I_y().get();
00787           send_info[d][num*3+2] = I.I_z().get();
00788           var_send[l][d][num] = Give_variable(I,l);
00789           num++;
00790         }
00791       }
00792       */
00793 
00794       // 5.2 part: send points
00795       num_message = 0;
00796       if(rank_source != -1 && num_rec_face_I[d]>0) {
00797         MPI_Irecv(receive_info[d] ,
00798                   3*num_rec_face_I[d],
00799                   MPI_INT,rank_source,40+d,comm,
00800                   &req[num_message]);
00801         ++num_message;
00802       }
00803       if(rank_destination != -1 && num_send_face_I[d]>0) {
00804         MPI_Isend(send_info[d],
00805                   3*num_send_face_I[d],
00806                   MPI_INT,rank_destination,40+d,comm,
00807                   &req[num_message]);
00808         ++num_message;
00809       }
00810       MPI_Waitall(num_message,req,status);
00811     
00812       // 5.3 part: set points
00813       for(num=0;num<num_rec_face_I[d];++num) {
00814         I = Index3D(receive_info[d][num*3],
00815                     receive_info[d][num*3+1],
00816                     receive_info[d][num*3+2]);
00817         var_receive[l][d][num] = Give_variable_slow(I,l);
00818         if(var_receive[l][d][num]==NULL) {
00819           Add_double(I,l,Give_max_num_var());
00820           var_receive[l][d][num] = Give_variable(I,l);
00821         }
00822 
00823 
00824 
00825         if(developer_version) {
00826           ptyp = Give_type(I);
00827           if(ptyp!=parallel_p && ptyp!=exterior) {
00828             cout << "\n error in Prepare_communication(), 5.3, my rank: " 
00829                  << my_rank 
00830                  << " Index: " 
00831                  << " x: " << receive_info[d][num*3]
00832                  << " y: " << receive_info[d][num*3+1]
00833                  << " z: " << receive_info[d][num*3+2]
00834                  << "typ: " << ptyp 
00835                  << endl;
00836           }
00837         }
00838         
00839 
00840 
00841 
00842            
00843 
00844            /* old
00845         var_receive[l][d][num] = 
00846           Give_variable(Index3D(receive_info[d][num*3],
00847                                 receive_info[d][num*3+1],
00848                                 receive_info[d][num*3+2]), 
00849                         l);
00850         if(developer_version) {
00851           if(var_receive[l][d][num]==NULL) {
00852             cout << " Error in 5.3. l: " << l 
00853                  << " Index: " 
00854                  << " x: " << receive_info[d][num*3]
00855                  << " y: " << receive_info[d][num*3+1]
00856                  << " z: " << receive_info[d][num*3+2]
00857                  << endl;
00858           }
00859         }
00860            */
00861       }
00862 
00863       if(l==Max_level()) {
00864         // 6.  Step: transmitter boundary points
00865         // 6.1 build informations about points to be send
00866         // How informations are stored:
00867         // Index of corner point, dir_3D       
00868         // ind_x, ind_y, ind_z     d      >   4 int
00869         num = 0;
00870         iterate_hash3 {
00871           if(Send_B_point_in_direction((dir_3D)d,bo2point,next_index)) {
00872             I = bo2point->ind;
00873             
00874             send_info[d][num*4]   = I.I_x().get();
00875             send_info[d][num*4+1] = I.I_y().get();
00876             send_info[d][num*4+2] = I.I_z().get();
00877             send_info[d][num*4+3] =  bo2point->direction;
00878             
00879             var_send[Max_level()][d][num_send_face_I[d]+num] = 
00880               Give_variable(I,(dir_3D)(bo2point->direction));
00881             num++;
00882           }
00883         }
00884     
00885         // 6.2 part: send points
00886         num_message = 0;
00887         if(rank_source != -1 && num_rec_face_B[d]>0) {
00888           MPI_Irecv(receive_info[d] ,
00889                     4*num_rec_face_B[d],
00890                     MPI_INT,rank_source,60+d,comm,
00891                     &req[num_message]);
00892           ++num_message;
00893         }
00894         if(rank_destination != -1 && num_send_face_B[d]>0) {
00895           MPI_Isend(send_info[d],
00896                     4*num_send_face_B[d],
00897                     MPI_INT,rank_destination,60+d,comm,
00898                     &req[num_message]);
00899           ++num_message;
00900         }
00901         MPI_Waitall(num_message,req,status);
00902         
00903         // 6.3 part: set points
00904         for(num=0;num<num_rec_face_B[d];++num) {
00905           if(developer_version) {
00906             if(Exists_Bo_Point(Index3D(receive_info[d][num*4],
00907                                        receive_info[d][num*4+1],
00908                                        receive_info[d][num*4+2]), 
00909                                (dir_3D)receive_info[d][num*4+3])==false) {
00910               cout << "\n error in Prepare_communication(), 6.3, \n my: " 
00911                    << my_rank << " my_I: "; 
00912               transform_coord(my_index.coordinate()).Print();
00913               cout << " d:"; Print((dir_3D)d); cout << endl;
00914               Index3D II;
00915               II = Index3D(receive_info[d][num*4],
00916                            receive_info[d][num*4+1],
00917                            receive_info[d][num*4+2]);
00918               II.Print();
00919               transform_coord(II.coordinate()).Print();
00920               cout << " direction: "; Print((dir_3D)receive_info[d][num*4+3]);
00921               cout << endl;
00922             }
00923           }
00924 
00925           var_receive[Max_level()][d][num_rec_face_I[d]+num] = 
00926             Give_variable(Index3D(receive_info[d][num*4],
00927                                   receive_info[d][num*4+1],
00928                                   receive_info[d][num*4+2]), 
00929                           (dir_3D)receive_info[d][num*4+3]);
00930         }
00931         
00932         // 7. Step: transmitter Zellfreiheitsgrade points
00933         // 7.1 build informations about points to be send
00934         // How informations are stored:
00935         // Index of point       
00936         // ind_x, ind_y, ind_z     >   3 int
00937         num = 0;
00938 
00939         iterate_hash2 {
00940           if(Send_Z_point_in_direction((dir_3D)d,bocell,next_index)) {
00941             I = bocell->ind;
00942             
00943             send_info[d][num*3]   = I.I_x().get();
00944             send_info[d][num*3+1] = I.I_y().get();
00945             send_info[d][num*3+2] = I.I_z().get();
00946             //  var_send[d][num] = Give_variable_cellpoi(I);
00947             var_send[Max_level()][d][num_send_face_I[d]+num_send_face_B[d]
00948                                     +num] = 
00949               bocell->var;
00950             num++;
00951           }
00952         }
00953         
00954         // 7.2 part: send points
00955         num_message = 0;
00956         if(rank_source != -1 && num_rec_face_Z[d]) {
00957           MPI_Irecv(receive_info[d] ,
00958                     3*num_rec_face_Z[d],
00959                     MPI_INT,rank_source,80+d,comm,
00960                     &req[num_message]);
00961           ++num_message;
00962         }
00963         if(rank_destination != -1 && num_send_face_Z[d]) {
00964           MPI_Isend(send_info[d],
00965                     3*num_send_face_Z[d],
00966                     MPI_INT,rank_destination,80+d,comm,
00967                     &req[num_message]);
00968           ++num_message;
00969         }
00970         MPI_Waitall(num_message,req,status);
00971       
00972         // 7.3 part: set points
00973         for(num=0;num<num_rec_face_Z[d];++num) {
00974           var_receive[Max_level()][d][num_rec_face_I[d]+num_rec_face_B[d]
00975                                      +num] = 
00976             Give_variable_cellpoi(Index3D(receive_info[d][num*3],
00977                                           receive_info[d][num*3+1],
00978                                           receive_info[d][num*3+2]));
00979           if(developer_version) {
00980             if(var_receive[Max_level()][d][num_rec_face_I[d]+num_rec_face_B[d]
00981                                           +num] == NULL) {
00982               cout << " Error in 7.3, Prepare_communication_all_grids()! "
00983                    << " my rank: " << my_rank
00984                    << " from: " << rank_source << endl;
00985             }
00986           }
00987         }
00988       }
00989       // 8. part delete 
00990       if(send_info[d]!=NULL) {
00991         delete(send_info[d]);
00992         send_info[d] = NULL;
00993       }
00994       if(receive_info[d]!=NULL) {
00995         delete(receive_info[d]);
00996         receive_info[d] = NULL;
00997       }
00998     }
00999     delete invariant_array;
01000   }
01001 }
01002   
01003 void Grid::Full_update_Variable(int level) {
01004   int* array_num;
01005   int max_var;
01006   
01007   max_var = Give_max_num_var();
01008   array_num  = new int[max_var];
01009   for(int i=0;i<max_var;++i) {
01010     array_num[i] = i;
01011   }
01012   Update_Variable(Wdir,array_num,max_var,level);
01013   Update_Variable(Edir,array_num,max_var,level);
01014   Update_Variable(Ndir,array_num,max_var,level);
01015   Update_Variable(Sdir,array_num,max_var,level);
01016   Update_Variable(Tdir,array_num,max_var,level);
01017   Update_Variable(Ddir,array_num,max_var,level);
01018   
01019   /*
01020   Index1D nextS;
01021   int d,i, num, lr, iv;
01022   MPI_Request req[12];
01023   MPI_Status status[12];
01024   int num_message;
01025   double*  send_info[6];
01026   double*  receive_info[6];
01027   int number_send[6];
01028   int number_receive[6];
01029 
01030   int number_vars;
01031 
01032   if(I_am_active()) {
01033     // 1. Step: calc number to be send and receive
01034     number_vars = Give_max_num_var();
01035     for(d=0;d<6;++d) {
01036       number_send[d] = 
01037         number_send_face[level][d] * number_vars;
01038       number_receive[d] = 
01039         number_receive_face[level][d] * number_vars;
01040     }
01041     
01042     // 2. Step: build informations about var(double) to be send  
01043     for(d=0;d<6;++d) {
01044       send_info[d]    = Give_buffer_send((dir_3D)d,number_send[d]); 
01045       receive_info[d] = Give_buffer_receive((dir_3D)d,number_receive[d]); 
01046     }
01047     
01048     for(i=0;i<3;++i) {
01049       for(lr=0;lr<2;++lr) {
01050         d =  Add((dir1D)lr,(main_dir_3D)i);
01051         for(num=0;num<number_send_face[level][d];++num) {
01052           for(iv=0;iv<Give_max_num_var();++iv) {
01053             send_info[d][num*number_vars+iv] = 
01054               var_send[level][d][num][iv];
01055             //    send_info[d][num*number_vars+iv] = 100.0;
01056           }
01057         }
01058       }
01059       
01060       // 4. Step: send doubles
01061       num_message = 0;
01062       for(lr=0;lr<2;++lr) {
01063         d =  Add((dir1D)lr,(main_dir_3D)i);
01064         if(give_next_rank_source((dir_3D)d) != -1 && number_receive[d]>0) {
01065           MPI_Irecv(receive_info[d] ,number_receive[d],
01066                     MPI_DOUBLE,give_next_rank_source((dir_3D)d),10,comm,
01067                     &req[num_message]);
01068           ++num_message;
01069         }
01070       }
01071       for(lr=0;lr<2;++lr) {
01072         d =  Add((dir1D)lr,(main_dir_3D)i);
01073         if(give_next_rank_destination((dir_3D)d) != -1 && number_send[d]>0) {
01074           MPI_Isend(send_info[d],number_send[d],
01075                     MPI_DOUBLE,give_next_rank_destination((dir_3D)d),10,comm,
01076                     &req[num_message]);
01077           ++num_message;
01078         }
01079       }
01080       MPI_Waitall(num_message,req,status);
01081       
01082       // 5. part: set informations about var(double)
01083       for(lr=0;lr<2;++lr) {
01084         d =  Add((dir1D)lr,(main_dir_3D)i);
01085         for(num=0;num<number_receive_face[level][d];++num) {
01086           for(iv=0;iv<Give_max_num_var();++iv) {
01087             var_receive[level][d][num][iv] = 
01088               receive_info[d][num*number_vars+iv];
01089             // var_receive[d][num][i] = 100;
01090           }
01091         }
01092       }
01093     }
01094   }
01095   */
01096 }
01097 
01098 void Grid::Update_Variable_no_full(int* array_num, int anz_num,
01099                                    int level) {
01100   for(int i=0;i<26;++i) {
01101     Update_Variable(i, array_num, anz_num, level);
01102   }
01103 }
01104 void Grid::Update_Variable(dir_3D d, int* array_num, int anz_num, int level) {
01105   if(d==Wdir) {
01106     Update_Variable((int)Wdir, array_num, anz_num, level);
01107 
01108     Update_Variable((int)(WDed+6), array_num, anz_num, level);
01109     Update_Variable((int)(WTed+6), array_num, anz_num, level);
01110     Update_Variable((int)(SWed+6), array_num, anz_num, level);
01111     Update_Variable((int)(NWed+6), array_num, anz_num, level);
01112 
01113     Update_Variable((int)(WSTd+18), array_num, anz_num, level);
01114     Update_Variable((int)(WSDd+18), array_num, anz_num, level);
01115     Update_Variable((int)(WSTd+18), array_num, anz_num, level);
01116     Update_Variable((int)(WNDd+18), array_num, anz_num, level);
01117     Update_Variable((int)(WNTd+18), array_num, anz_num, level);
01118   }
01119   if(d==Edir) {
01120     Update_Variable((int)Edir, array_num, anz_num, level);
01121 
01122     Update_Variable((int)(EDed+6), array_num, anz_num, level);
01123     Update_Variable((int)(ETed+6), array_num, anz_num, level);
01124     Update_Variable((int)(SEed+6), array_num, anz_num, level);
01125     Update_Variable((int)(NEed+6), array_num, anz_num, level);
01126 
01127     Update_Variable((int)(ESDd+18), array_num, anz_num, level);
01128     Update_Variable((int)(ESTd+18), array_num, anz_num, level);
01129     Update_Variable((int)(ENDd+18), array_num, anz_num, level);
01130     Update_Variable((int)(ENTd+18), array_num, anz_num, level);
01131   }
01132   if(d==Ndir) {
01133     Update_Variable((int)Ndir, array_num, anz_num, level);
01134 
01135     Update_Variable((int)(NDed+6), array_num, anz_num, level);
01136     Update_Variable((int)(NTed+6), array_num, anz_num, level);
01137     Update_Variable((int)(NWed+6), array_num, anz_num, level);
01138     Update_Variable((int)(NEed+6), array_num, anz_num, level);
01139 
01140     Update_Variable((int)(ENTd+18), array_num, anz_num, level);
01141     Update_Variable((int)(WNTd+18), array_num, anz_num, level);
01142     Update_Variable((int)(ENDd+18), array_num, anz_num, level);
01143     Update_Variable((int)(WNDd+18), array_num, anz_num, level);
01144   }
01145   if(d==Sdir) {
01146     Update_Variable((int)Sdir, array_num, anz_num, level);
01147 
01148     Update_Variable((int)(SDed+6), array_num, anz_num, level);
01149     Update_Variable((int)(STed+6), array_num, anz_num, level);
01150     Update_Variable((int)(SWed+6), array_num, anz_num, level);
01151     Update_Variable((int)(SEed+6), array_num, anz_num, level);
01152 
01153     Update_Variable((int)(ESTd+18), array_num, anz_num, level);
01154     Update_Variable((int)(WSTd+18), array_num, anz_num, level);
01155     Update_Variable((int)(ESDd+18), array_num, anz_num, level);
01156     Update_Variable((int)(WSDd+18), array_num, anz_num, level);
01157   }
01158   if(d==Ddir) {
01159     Update_Variable((int)Ddir, array_num, anz_num, level);
01160 
01161     Update_Variable((int)(SDed+6), array_num, anz_num, level);
01162     Update_Variable((int)(NDed+6), array_num, anz_num, level);
01163     Update_Variable((int)(WDed+6), array_num, anz_num, level);
01164     Update_Variable((int)(EDed+6), array_num, anz_num, level);
01165 
01166     Update_Variable((int)(ESDd+18), array_num, anz_num, level);
01167     Update_Variable((int)(WSDd+18), array_num, anz_num, level);
01168     Update_Variable((int)(ENDd+18), array_num, anz_num, level);
01169     Update_Variable((int)(WNDd+18), array_num, anz_num, level);
01170   }
01171   if(d==Tdir) {
01172     Update_Variable((int)Tdir, array_num, anz_num, level);
01173 
01174     Update_Variable((int)(STed+6), array_num, anz_num, level);
01175     Update_Variable((int)(NTed+6), array_num, anz_num, level);
01176     Update_Variable((int)(WTed+6), array_num, anz_num, level);
01177     Update_Variable((int)(ETed+6), array_num, anz_num, level);
01178 
01179     Update_Variable((int)(ESTd+18), array_num, anz_num, level);
01180     Update_Variable((int)(WSTd+18), array_num, anz_num, level);
01181     Update_Variable((int)(ENTd+18), array_num, anz_num, level);
01182     Update_Variable((int)(WNTd+18), array_num, anz_num, level);
01183   }
01184 }
01185 
01186 void Grid::Update_Variable(int d, int* array_num, int anz_num, int level) {
01187   Index1D nextS;
01188   int num, iv;
01189   MPI_Request req[12];
01190   MPI_Status status[12];
01191   int num_message;
01192   double*  send_info;
01193   double*  receive_info;
01194   int number_send;
01195   int number_receive;
01196 
01197   //  int number_vars;
01198   double** var_sendd;
01199   double** var_received;
01200 
01201   int rank_source, rank_destination;
01202 
01203   /*
01204   MPI_Barrier(MPI_COMM_WORLD);
01205   if(my_rank==0) {
01206     cout << " Parallelisierung in Richtung: " << d << endl;
01207     for(iv=0;iv<anz_num;++iv) {
01208       cout << " Nummer: " << array_num[iv] << endl;
01209     }
01210   }
01211   MPI_Barrier(MPI_COMM_WORLD);
01212   */
01213 
01214   //if(I_am_active() && Max_level()==level) {
01215   if(I_am_active() && anz_num != 0) {
01216     // 0. part: basic infos about send and recieve
01217    if(d<6) {
01218       rank_source      = give_next_rank_source((dir_3D)d,level);
01219       rank_destination = give_next_rank_destination((dir_3D)d,level);  
01220     }
01221     else if(d<18) {
01222       rank_source      = give_next_rank_source((Edges_cell)(d-6),level);
01223       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level);  
01224     }
01225     else {
01226       rank_source      = give_next_rank_source((dir_sons)(d-18),level);
01227       rank_destination = give_next_rank_destination((dir_sons)(d-18),level);
01228     }
01229 
01230     // 1. Step: calc number to be send and receive
01231     number_send = 
01232       number_send_face[level][d] * anz_num;
01233     number_receive = 
01234       number_receive_face[level][d] * anz_num;
01235    
01236     // 2. Step: build informations about var(double) to be send  
01237     send_info    = Give_buffer_send(d,number_send); 
01238     receive_info = Give_buffer_receive(d,number_receive); 
01239 
01240     var_sendd=var_send[level][d];
01241     for(num=0;num<number_send_face[level][d];++num) {
01242       for(iv=0;iv<anz_num;++iv) {
01243         send_info[num*anz_num+iv] = 
01244           var_sendd[num][array_num[iv]];
01245       }
01246     }
01247 
01248     // 4. Step: send doubles
01249     num_message = 0;
01250     if(rank_source != -1 && number_receive>0) {
01251       MPI_Irecv(receive_info ,number_receive,
01252                 MPI_DOUBLE,rank_source,26,comm,
01253                 &req[num_message]);
01254       ++num_message;
01255     }      
01256     if(rank_destination != -1 && number_send>0) {
01257       MPI_Isend(send_info,number_send,
01258                 MPI_DOUBLE,rank_destination,26,comm,
01259                 &req[num_message]);
01260       ++num_message;
01261     }
01262     MPI_Waitall(num_message,req,status);
01263 
01264     // 5. part: set informations about var(double)
01265     var_received=var_receive[level][d];
01266     for(num=0;num<number_receive_face[level][d];++num) {
01267       for(iv=0;iv<anz_num;++iv) {
01268         var_received[num][array_num[iv]] = 
01269           receive_info[num*anz_num+iv];
01270       }
01271     }
01272   }
01273 }
01274 
01275 void Grid::Update_Variable_for_prolongation(int* array_num, int anz_num,
01276                                             int level) {
01277   for(int i=0;i<26;++i) {
01278     Update_Variable_for_prolongation(i, array_num, anz_num, level);
01279   }
01280 }
01281 
01282 void Grid::Update_Variable_for_prolongation(dir_3D d, int* array_num, int anz_num, int level) {
01283   if(d==Wdir) {
01284     Update_Variable_for_prolongation((int)Wdir, array_num, anz_num, level);
01285 
01286     Update_Variable_for_prolongation((int)(WDed+6), array_num, anz_num, level);
01287     Update_Variable_for_prolongation((int)(WTed+6), array_num, anz_num, level);
01288     Update_Variable_for_prolongation((int)(SWed+6), array_num, anz_num, level);
01289     Update_Variable_for_prolongation((int)(NWed+6), array_num, anz_num, level);
01290 
01291     Update_Variable_for_prolongation((int)(WSTd+18), array_num, anz_num, level);
01292     Update_Variable_for_prolongation((int)(WSDd+18), array_num, anz_num, level);
01293     Update_Variable_for_prolongation((int)(WSTd+18), array_num, anz_num, level);
01294     Update_Variable_for_prolongation((int)(WNDd+18), array_num, anz_num, level);
01295     Update_Variable_for_prolongation((int)(WNTd+18), array_num, anz_num, level);
01296   }
01297   if(d==Edir) {
01298     Update_Variable_for_prolongation((int)Edir, array_num, anz_num, level);
01299 
01300     Update_Variable_for_prolongation((int)(EDed+6), array_num, anz_num, level);
01301     Update_Variable_for_prolongation((int)(ETed+6), array_num, anz_num, level);
01302     Update_Variable_for_prolongation((int)(SEed+6), array_num, anz_num, level);
01303     Update_Variable_for_prolongation((int)(NEed+6), array_num, anz_num, level);
01304 
01305     Update_Variable_for_prolongation((int)(ESDd+18), array_num, anz_num, level);
01306     Update_Variable_for_prolongation((int)(ESTd+18), array_num, anz_num, level);
01307     Update_Variable_for_prolongation((int)(ENDd+18), array_num, anz_num, level);
01308     Update_Variable_for_prolongation((int)(ENTd+18), array_num, anz_num, level);
01309   }
01310   if(d==Ndir) {
01311     Update_Variable_for_prolongation((int)Ndir, array_num, anz_num, level);
01312 
01313     Update_Variable_for_prolongation((int)(NDed+6), array_num, anz_num, level);
01314     Update_Variable_for_prolongation((int)(NTed+6), array_num, anz_num, level);
01315     Update_Variable_for_prolongation((int)(NWed+6), array_num, anz_num, level);
01316     Update_Variable_for_prolongation((int)(NEed+6), array_num, anz_num, level);
01317 
01318     Update_Variable_for_prolongation((int)(ENTd+18), array_num, anz_num, level);
01319     Update_Variable_for_prolongation((int)(WNTd+18), array_num, anz_num, level);
01320     Update_Variable_for_prolongation((int)(ENDd+18), array_num, anz_num, level);
01321     Update_Variable_for_prolongation((int)(WNDd+18), array_num, anz_num, level);
01322   }
01323   if(d==Sdir) {
01324     Update_Variable_for_prolongation((int)Sdir, array_num, anz_num, level);
01325 
01326     Update_Variable_for_prolongation((int)(SDed+6), array_num, anz_num, level);
01327     Update_Variable_for_prolongation((int)(STed+6), array_num, anz_num, level);
01328     Update_Variable_for_prolongation((int)(SWed+6), array_num, anz_num, level);
01329     Update_Variable_for_prolongation((int)(SEed+6), array_num, anz_num, level);
01330 
01331     Update_Variable_for_prolongation((int)(ESTd+18), array_num, anz_num, level);
01332     Update_Variable_for_prolongation((int)(WSTd+18), array_num, anz_num, level);
01333     Update_Variable_for_prolongation((int)(ESDd+18), array_num, anz_num, level);
01334     Update_Variable_for_prolongation((int)(WSDd+18), array_num, anz_num, level);
01335   }
01336   if(d==Ddir) {
01337     Update_Variable_for_prolongation((int)Ddir, array_num, anz_num, level);
01338 
01339     Update_Variable_for_prolongation((int)(SDed+6), array_num, anz_num, level);
01340     Update_Variable_for_prolongation((int)(NDed+6), array_num, anz_num, level);
01341     Update_Variable_for_prolongation((int)(WDed+6), array_num, anz_num, level);
01342     Update_Variable_for_prolongation((int)(EDed+6), array_num, anz_num, level);
01343 
01344     Update_Variable_for_prolongation((int)(ESDd+18), array_num, anz_num, level);
01345     Update_Variable_for_prolongation((int)(WSDd+18), array_num, anz_num, level);
01346     Update_Variable_for_prolongation((int)(ENDd+18), array_num, anz_num, level);
01347     Update_Variable_for_prolongation((int)(WNDd+18), array_num, anz_num, level);
01348   }
01349   if(d==Tdir) {
01350     Update_Variable_for_prolongation((int)Tdir, array_num, anz_num, level);
01351 
01352     Update_Variable_for_prolongation((int)(STed+6), array_num, anz_num, level);
01353     Update_Variable_for_prolongation((int)(NTed+6), array_num, anz_num, level);
01354     Update_Variable_for_prolongation((int)(WTed+6), array_num, anz_num, level);
01355     Update_Variable_for_prolongation((int)(ETed+6), array_num, anz_num, level);
01356 
01357     Update_Variable_for_prolongation((int)(ESTd+18), array_num, anz_num, level);
01358     Update_Variable_for_prolongation((int)(WSTd+18), array_num, anz_num, level);
01359     Update_Variable_for_prolongation((int)(ENTd+18), array_num, anz_num, level);
01360     Update_Variable_for_prolongation((int)(WNTd+18), array_num, anz_num, level);
01361   }
01362 }
01363 
01364 void Grid::Update_Variable_for_prolongation(int d, 
01365                                             int* array_num, 
01366                                             int anz_num, int level) {
01367   Index1D nextS;
01368   int num, iv;
01369   MPI_Request req[12];
01370   MPI_Status status[12];
01371   int num_message;
01372   double*  send_info;
01373   double*  receive_info;
01374   int number_send;
01375   int number_receive;
01376 
01377   //  int number_vars;
01378   double** var_sendd;
01379   double** var_received;
01380 
01381   int rank_source, rank_destination;
01382 
01383   // level is level with respect to coarse points
01384   // send point level is 1 larger
01385   if(I_am_active()) {
01386     // 0. part: basic infos about send and recieve
01387     if(d<6) {
01388       rank_source      = give_next_rank_source((dir_3D)d,level+1);
01389       rank_destination = give_next_rank_destination((dir_3D)d,level+1);  
01390     }
01391     else if(d<18) {
01392       rank_source      = give_next_rank_source((Edges_cell)(d-6),level+1);
01393       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level+1);  
01394     }
01395     else {
01396       rank_source      = give_next_rank_source((dir_sons)(d-18),level+1);
01397       rank_destination = give_next_rank_destination((dir_sons)(d-18),level+1);
01398     }
01399 
01400     // 1. Step: calc number to be send and receive
01401     number_send = 
01402       number_send_face_prolongation[level][d] * anz_num;
01403     number_receive = 
01404       number_receive_face_prolongation[level][d] * anz_num;
01405    
01406     // 2. Step: build informations about var(double) to be send  
01407     send_info    = Give_buffer_send(d,number_send); 
01408     receive_info = Give_buffer_receive(d,number_receive); 
01409 
01410     var_sendd=var_send_prolongation[level][d];
01411     for(num=0;num<number_send_face_prolongation[level][d];++num) {
01412       for(iv=0;iv<anz_num;++iv) {
01413         send_info[num*anz_num+iv] = 
01414           var_sendd[num][array_num[iv]];
01415       }
01416     }
01417 
01418     // 4. Step: send doubles
01419     num_message = 0;
01420     if(rank_source != -1 && number_receive>0) {
01421       MPI_Irecv(receive_info ,number_receive,
01422                 MPI_DOUBLE,rank_source,130+d,comm,
01423                 &req[num_message]);
01424       ++num_message;
01425     }      
01426     if(rank_destination != -1 && number_send>0) {
01427       MPI_Isend(send_info,number_send,
01428                 MPI_DOUBLE,rank_destination,130+d,comm,
01429                 &req[num_message]);
01430       ++num_message;
01431     }
01432     MPI_Waitall(num_message,req,status);
01433 
01434     // 5. part: set informations about var(double)
01435     var_received=var_receive_prolongation[level][d];
01436     for(num=0;num<number_receive_face_prolongation[level][d];++num) {
01437       for(iv=0;iv<anz_num;++iv) {
01438         var_received[num][array_num[iv]] = 
01439           receive_info[num*anz_num+iv];
01440       }
01441     }
01442   }
01443 }
01444 
01445 
01446 
01447 void Grid::Print_number_of_communication_doubles() {
01448   int sum, total_sum;
01449   int rank_destination, d;
01450   if(my_rank==0)
01451     cout << "\n Print number of communication doubles: " 
01452          << endl;
01453   for(int l=min_level;l<=max_level;++l) {
01454     sum = 0;
01455     if(I_am_active()) {
01456       for(d=0;d<26;++d) {
01457         if(d<6) {
01458           rank_destination = give_next_rank_destination((dir_3D)d,l);  
01459         }
01460         else if(d<18) {
01461           rank_destination = give_next_rank_destination((Edges_cell)(d-6),l);
01462         }
01463         else {
01464           rank_destination = give_next_rank_destination((dir_sons)(d-18),l);
01465         }
01466         
01467         if(rank_destination != -1)
01468           sum = sum + number_send_face[l][d];
01469       }
01470     }
01471     MPI_Allreduce(&sum,&total_sum,1,MPI_INT,MPI_SUM,comm);
01472     if(my_rank==0)
01473       cout << " level: " << l << " is: " << total_sum << endl;
01474   }
01475 }
01476 
01477 
01478 
01479 
01480 
01481 
01482 
01484 
01485 void Grid::Sum_ghost_nodes_Variable(int* array_num, int anz_num, int level) {
01486   for(int d=0;d<26;++d)
01487     Sum_ghost_nodes_Variable(d, array_num, anz_num, level);
01488 }
01489 
01490 void Grid::Sum_ghost_nodes_Variable(int d, 
01491                                     int* array_num, int anz_num, int level) {
01492   Index1D nextS;
01493   int num, iv;
01494   MPI_Request req[12];
01495   MPI_Status status[12];
01496   int num_message;
01497   double*  send_info;
01498   double*  receive_info;
01499   int number_send;
01500   int number_receive;
01501 
01502   //  int number_vars;
01503   double** var_sendd;
01504   double** var_received;
01505 
01506   int rank_source, rank_destination;
01507 
01508 
01509   //if(I_am_active() && Max_level()==level) {
01510   if(I_am_active() && anz_num != 0) {
01511     // 0. part: basic infos about send and recieve
01512    if(d<6) {
01513       rank_source      = give_next_rank_source((dir_3D)d,level);
01514       rank_destination = give_next_rank_destination((dir_3D)d,level);  
01515     }
01516     else if(d<18) {
01517       rank_source      = give_next_rank_source((Edges_cell)(d-6),level);
01518       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level);  
01519     }
01520     else {
01521       rank_source      = give_next_rank_source((dir_sons)(d-18),level);
01522       rank_destination = give_next_rank_destination((dir_sons)(d-18),level);
01523     }
01524 
01525     // 1. Step: calc number to be send and receive
01526     number_send = 
01527       number_send_face[level][d] * anz_num;
01528     number_receive = 
01529       number_receive_face[level][d] * anz_num;
01530    
01531     // 2. Step: build informations about var(double) to be send  
01532     send_info    = Give_buffer_send(d,number_send); 
01533     receive_info = Give_buffer_receive(d,number_receive); 
01534 
01535     var_received=var_receive[level][d];
01536     for(num=0;num<number_receive_face[level][d];++num) {
01537       for(iv=0;iv<anz_num;++iv) {
01538         receive_info[num*anz_num+iv] =
01539         var_received[num][array_num[iv]];
01540       }
01541     }
01542 
01543     // 4. Step: send doubles
01544     num_message = 0;
01545     if(rank_destination != -1 && number_send>0) {
01546       MPI_Irecv(send_info,number_send,
01547                 MPI_DOUBLE,rank_destination,26,comm,
01548                 &req[num_message]);
01549       ++num_message;
01550     }
01551     if(rank_source != -1 && number_receive>0) {
01552       MPI_Isend(receive_info ,number_receive,
01553                 MPI_DOUBLE,rank_source,26,comm,
01554                 &req[num_message]);
01555       ++num_message;
01556     }      
01557 
01558     MPI_Waitall(num_message,req,status);
01559 
01560     // 5. part: set informations about var(double)
01561     if(rank_destination != -1 && number_send>0) {
01562       var_sendd=var_send[level][d];
01563       for(num=0;num<number_send_face[level][d];++num) {
01564         for(iv=0;iv<anz_num;++iv) {
01565           var_sendd[num][array_num[iv]] +=
01566             send_info[num*anz_num+iv];
01567         }
01568       } 
01569     }
01570   }
01571 }

Generated on Fri Nov 2 01:25:56 2007 for IPPL by doxygen 1.3.5