src/expde/grid/labbo.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 // labbo.cc
00020 //
00021 // ------------------------------------------------------------
00022 
00023 
00024 // Include level 0:
00025 #ifdef COMP_GNUOLD
00026  #include <iostream.h>
00027  #include <fstream.h>
00028  #include <time.h>
00029  #include <math.h>
00030 #else
00031  #include <iostream>
00032  #include <fstream> 
00033  #include <ctime>
00034  #include <cmath>
00035 #endif 
00036 
00037 
00038 // Include level 0:
00039 #include "../parser.h"
00040 
00041 // Include level 1:
00042 #include "../paramete.h"
00043 #include "../abbrevi.h"
00044 #include "../math_lib/math_lib.h"
00045 
00046 // Include level 2:
00047 #include "../basic/basic.h"
00048 
00049 // Include level 3:
00050 #include "../domain/domain.h"
00051 
00052 // Include level 4:
00053 #include "../formulas/boundy.h"
00054 #include "../formulas/loc_sten.h"
00055 
00056 // Include level 5:
00057 #include "gpar.h"
00058 #include "parallel.h"
00059 #include "mgcoeff.h"
00060 #include "sto_man.h"
00061 #include "gridbase.h"
00062 #include "grid.h"
00063 
00065 // Memberfunctions for boundary label
00066 // 1. at boundary points
00067 // 2. for multigrid points
00068 // 3. Parallelization
00070 
00071 bool jetzt;
00072 
00073 
00074 // 1. at boundary points
00076 
00077 
00078 int Storage_manager::Give_number_label_bo() {
00079   int i;
00080   for(i=0;i<Max_label_number;++i) {
00081     if(occupied[i]==false) {
00082       occupied[i]=true;
00083       if(i>=max_num_label) {
00084         // erhoehe max_num_label
00085         increase_max_num_label = true;
00086         ++max_num_label;
00087       }
00088       else 
00089         // erhoehe nicht max_num_label
00090         increase_max_num_label = false;
00091       return i+1;
00092     }
00093   }
00094   cout << " Error: Maximal number of labels for the boundary reached! "
00095        << endl;
00096   return 0;
00097 }
00098 
00099 void Storage_manager::Remove_number_label_bo(int num) {
00100   if(developer_version) {
00101     if(num>Max_label_number || num<=0) {
00102       cout << " Label number wrong in Grid_base::Remove_number_label_bo! "
00103            << endl;
00104       return;
00105     }
00106   }
00107   occupied[num-1]=false;
00108 }
00109 
00110 bool Grid_base::Give_label_bo(Index3D Ind, dir_3D dir, int num) 
00111      const {
00112   static Point_hashtable3* point;
00113   point = hashtable3_start
00114     [hashtable3_function(Ind.ind_x.index,Ind.ind_y.index,
00115                          Ind.ind_z.index, dir, hashtable3_leng)];
00116   while(point!=NULL) {
00117     if(point->ind==Ind && point->direction==dir) return point->Give_label(num);
00118     point = point->next;
00119   }
00120   return point->Give_label(num);
00121 }
00122 
00123 
00124 /*
00125 void Grid_base::Put_label_bo(bool lab, 
00126                              Index3D Ind, dir_3D dir, int num) 
00127      const {
00128   static Point_hashtable3* point;
00129   point = hashtable3_start
00130     [hashtable3_function(Ind.ind_x.index,Ind.ind_y.index,
00131                          Ind.ind_z.index, dir, hashtable3_leng)];
00132   while(point!=NULL) {
00133     if(point->ind==Ind && point->direction==dir) point->Put_label(lab,num);
00134     point = point->next;
00135   }
00136   point->Put_label(lab,num);
00137 }
00138 */
00139 
00140 void Grid_base::Put_label_bo(bool (*Formula)(double x,double y,double z),
00141                             int num) {
00142   D3vector Vl, Vm, Vr;
00143   iterate_hash3 {
00144     Vl=transform_coord_org_shift(bo2point->ind,bo2point->direction,-shift_eps);
00145     Vm=transform_coord_org_shift(bo2point->ind,bo2point->direction,0);
00146     Vr=transform_coord_org_shift(bo2point->ind,bo2point->direction,shift_eps);
00147     bo2point->Put_label(Formula(Vl.x,Vl.y,Vl.z)
00148                         || Formula(Vm.x,Vm.y,Vm.z)
00149                         || Formula(Vr.x,Vr.y,Vr.z),num);
00150   }
00151   Restrict_label_bo(num);
00152 }
00153 
00154 void Grid_base::Put_label_bo(int num_var, int num) {
00155   double epsi;
00156   epsi = finest_mesh_size() * eps_boundary_var;
00157   iterate_hash3 {
00158     bo2point->Put_label( bo2point->var[num_var] > epsi, num);
00159   }
00160   Restrict_label_bo(num);
00161 }
00162 
00163 void Grid_base::Put_label_bo_complementary(int num_org, int num) {
00164   iterate_hash3 {
00165     if(bo2point->Give_label(num_org)==true) 
00166       bo2point->Put_label(false,num);
00167     else 
00168       bo2point->Put_label(true,num);
00169   }
00170   Restrict_label_bo(num);
00171 }
00172 
00173 void Grid_base::Put_union_label(int num, int numa, int numb) {
00174   iterate_hash3 {
00175     bo2point->Put_label(bo2point->Give_label(numa) ||
00176                         bo2point->Give_label(numb),num);
00177   }
00178   Restrict_label_bo(num);
00179 }
00180 
00181 
00182 // 2. for multigrid points
00184 
00185 /*
00186   How to construct and restrict label_bo:
00187   1. label_bo must be set at all corners of boundary cells
00188      on all levels, since res_op.h needs these informations
00189   2. on finest grid:
00190           - boundary points -> true or false
00191           - set all fine_bocells -> false
00192           - Zellfreiheitgrade, if one corner Dirichlet, then all of them
00193   3. recursion: all corners of edges, faces and volumnes
00194                 must be set.
00195                 -> true if Dirichlet and used for interpolation
00196                 -> nothing else
00197      by this all coarse bocells have a label_bo
00198  */
00199 
00200 void Grid_base::Restrict_label_bo(int num) {
00201   Index3D I, Inext;
00202   int i;
00203   bool Dirichlet;
00204   Pointtype type;
00205 
00206   if(I_am_active()) {
00207     // 1. Schritt: increase labels of multigrid points
00208     //             if max_num_label>1 and increase_max_num_label
00209     int *new_l, *old_l;
00210     // for multigrid points
00211     if(Give_max_num_label()>1 && Give_increase_max_num_label()) {
00212       iterate_hash1 {
00213         I = point1->Give_Index();
00214         type = Give_type(I);
00215         if(type==multigrid || type==parallel_p) {
00216           old_l = point1->label_level;
00217           if(old_l!=NULL) {
00218             new_l = new int[Give_max_num_label()];
00219             for(i=0;i<Give_max_num_label()-1;++i)
00220               new_l[i]=old_l[i];
00221             new_l[Give_max_num_label()-1] = 0;
00222             delete old_l;
00223             point1->label_level = new_l;
00224           }
00225         }
00226       }
00227     }
00228     
00229     if(Give_increase_max_num_label()) {
00230       // for Zellfreiheitgrade
00231       iterate_hash2 {
00232         if(bocell->Exists_bocellpoint()) {
00233           old_l = bocell->Give_label_level();
00234           new_l = new int[Give_max_num_label()];
00235           for(i=0;i<Give_max_num_label()-1;++i)
00236             new_l[i]=old_l[i];
00237           new_l[Give_max_num_label()-1] = 0;
00238           bocell->Set_label_level(new_l);
00239         }
00240       }
00241     }
00242     
00243     // 2. Schritt: study points on finest grid
00244     
00245     
00246     // for Rand first: dont_send
00247     iterate_hash1 {
00248       point1->send_label = dont_send;
00249     }
00250     
00251     
00252     // boundary points
00253     iterate_hash3 {
00254       I = bo2point->ind.next(bo2point->direction,max_level);
00255       // means: since one cannot run that boundary point -> Dirichlet
00256       Dirichlet = bo2point->Give_label(num)==false;
00257       type = Give_type(I);
00258       if(type==multigrid || type==parallel_p) {
00259         Set_label_bo_mg(I,num,Dirichlet,max_level);
00260       }
00261       Interpolation_label_bo(I,max_level,num,Dirichlet);
00262       Interpolation_label_bo(bo2point->ind,max_level,num,Dirichlet);
00263     }
00264     
00265     /*
00266       Send_label_parallel(max_level,num);
00267       Send_label_parallel(max_level-1,num);
00268       // for Rand first: dont_send
00269       iterate_hash1 {
00270       point1->send_label = dont_send;
00271       }
00272     */
00273     
00274     // Zellfreiheitgrade and all fine_bo_cells
00275     iterate_hash2 {
00276       I = bocell->Give_Index();
00277       // for all fine_bo_cells all corners must have a label_bo
00278       for(i=0;i<8;++i) {
00279         Inext = I.neighbour((dir_sons)i);
00280         type = Give_type(Inext);
00281         if(type==multigrid || type==parallel_p) {
00282         Set_label_bo_mg(Inext,num,false,max_level);
00283         }
00284       }
00285     }
00286     
00287     // Send_label_parallel(max_level,num);
00288     if(parallel_version)
00289       Send_label_parallel(max_level-1,num);
00290     
00291     /*
00292     // for Rand first: dont_send
00293     iterate_hash1 {
00294       point1->send_label = dont_send;
00295     }
00296     */
00297 
00298     iterate_hash2 {
00299       // Zellfreiheitgrad
00300       if(bocell->Exists_bocellpoint()) {
00301         Dirichlet = false;
00302         // calc Dirichlet cell or not
00303         for(i=0;i<8;++i) {
00304           Inext = I.neighbour((dir_sons)i);
00305           type = Give_type(Inext);
00306           if(type==multigrid || type==parallel_p)
00307             if(Give_label_bo_D(Inext,num,max_level)==true)
00308               Dirichlet = true;
00309         }
00310         // set type D-N at neighbour points
00311         for(i=0;i<8;++i) {
00312           Inext = I.neighbour((dir_sons)i);
00313           type = Give_type(Inext);
00314           if(type==multigrid || type==parallel_p) {
00315             Set_label_bo_mg(Inext,num,Dirichlet,max_level);
00316             Interpolation_label_bo(Inext, max_level, num, Dirichlet);
00317           }
00318         }
00319         if(Dirichlet)
00320           bocell->Set_Dirichlet(num);
00321       }
00322     }
00323     
00324     if(parallel_version)
00325       Send_label_parallel(max_level-1,num);
00326     
00327     // 3. Schritt: rekusive from fine to coarse
00328     int level;
00329     Index3D my_lev_index;
00330     
00331     for(level = max_level-1;
00332         level >= Give_my_coarsest_level() && level >1;
00333         --level) {
00334 
00335       /*
00336       // for Rand first: dont_send
00337       iterate_hash1 {
00338         point1->send_label = dont_send;
00339       }
00340       */
00341 
00342       my_lev_index = Give_my_level_index(level+1);
00343       iterate_hash1 {
00344         I = point1->Give_Index();
00345         if(I << my_lev_index)
00346           //if(point1->Give_typ()==multigrid && I.Tiefe()==level) {
00347           if(point1->global_typ==multigrid && I.Tiefe()==level) {
00348             if(point1->finest_level>=level) {
00349               if(point1->label_level!=NULL) {
00350                 Dirichlet = point1->Give_label(num)>=level;
00351                 Interpolation_label_bo(I, level, num, Dirichlet);
00352               }
00353               else {
00354                 cout << " error 1 in labbo.cc " << endl;
00355               }
00356             }
00357           }
00358       }
00359       
00360      
00361       if(parallel_version) {
00362         //      Send_label_parallel(level,num);
00363         Send_label_parallel(level-1,num);
00364       }
00365     }
00366     
00367     // 5. Schritt: 
00368     No_increase_max_num_label();
00369 
00370     // 6.Schritt:  Test ob alles ok.
00371     iterate_hash1 {
00372       if(point1->global_typ==multigrid &&
00373          point1->ind.Tiefe()<max_level && 
00374          point1->label_level==NULL) {
00375         cout << " error 2 in labbo.cc "
00376                  << " x: " << point1->ind.I_x().get()
00377                  << " y: " << point1->ind.I_y().get()
00378                  << " z: " << point1->ind.I_z().get()
00379                  << " my_rank: " << my_rank
00380                  << endl;
00381         point1->global_typ = exterior;
00382       }
00383     }
00384   }
00385 }
00386 
00387 
00388 void Grid_base::Interpolation_label_bo(Index3D I, int level, int num, 
00389                                        bool Dirichlet) {
00390 
00391 
00392 
00393 /*
00394 if((my_rank==9) && (I == Index3D(0,4,4))) {
00395   cout <<  " \n \n jezt \n \n " << endl;
00396 }
00397 
00398 if((my_rank==15) && (I == Index3D(0,4,8))) {
00399   cout <<  " \n \n jezt2 \n \n " << endl;
00400 }
00401 if((my_rank==15) && (I == Index3D(64,4,67))) {
00402   cout <<  " \n \n jezt3 \n \n " << endl;
00403 }
00404 */
00405 
00406   if(I.I_x().Tiefe()==level) {
00407     if(I.I_y().Tiefe()==level) {
00408       if(I.I_z().Tiefe()==level) {
00409         // middle point
00410         // edge between EST-WND
00411         Set_label_bo_mg(I.next_EST(level),num,Dirichlet,level-1);
00412         Set_label_bo_mg(I.next_WND(level),num,Dirichlet,level-1);
00413       }
00414       else {
00415         // TD-surface
00416         // edge between ES-WN
00417         Set_label_bo_mg(I.next_ES(level),num,Dirichlet,level-1);
00418         Set_label_bo_mg(I.next_WN(level),num,Dirichlet,level-1);
00419       }
00420     }
00421     else {
00422       if(I.I_z().Tiefe()==level) {
00423         // NS-surface
00424         // edge between WT-ED
00425         Set_label_bo_mg(I.next_WT(level),num,Dirichlet,level-1);
00426         Set_label_bo_mg(I.next_ED(level),num,Dirichlet,level-1);
00427       }
00428       else {
00429         // EW-edge
00430         Set_label_bo_mg(I.next_E(level),num,Dirichlet,level-1);
00431         Set_label_bo_mg(I.next_W(level),num,Dirichlet,level-1);
00432       }
00433     }
00434   }
00435   else {
00436     if(I.I_y().Tiefe()==level) {
00437       if(I.I_z().Tiefe()==level) {
00438         // EW-surface
00439         // edge between ST-ND
00440         Set_label_bo_mg(I.next_ST(level),num,Dirichlet,level-1);
00441         Set_label_bo_mg(I.next_ND(level),num,Dirichlet,level-1);
00442       }
00443       else {
00444         // NS-edge
00445         Set_label_bo_mg(I.next_N(level),num,Dirichlet,level-1);
00446         Set_label_bo_mg(I.next_S(level),num,Dirichlet,level-1);
00447       }
00448     }
00449     else {
00450       if(I.I_z().Tiefe()==level) {
00451         // TD-edge
00452         Set_label_bo_mg(I.next_T(level),num,Dirichlet,level-1);
00453         Set_label_bo_mg(I.next_D(level),num,Dirichlet,level-1);
00454       }
00455     }
00456   }
00457 }
00458 
00459 
00460 
00461 void Grid_base::Set_label_bo_mg(Index3D I, int num, bool label, int level) {
00462   Point_hashtable1* point;
00463   Pointtype type;
00464 
00465 
00466   point=hashtable1_point(I);
00467 
00468   if(developer_version) {
00469     if(point==NULL) {
00470       cout << " Error 1 in Set_label_bo_mg(...)"
00471            << " my_rank: " << my_rank 
00472            << " x: " << I.I_x().get()
00473            << " y: " << I.I_y().get()
00474            << " z: " << I.I_z().get()
00475            << endl;
00476     }
00477   }
00478 
00479   /*
00480 if(I==Index3D(9,8,16)) {
00481   //if(I==Index3D(4,17,16) && my_rank==17) {
00482   cout << " \n \n set: level: " 
00483        << level 
00484        << " my rank: " << my_rank << endl;
00485 
00486   // point=NULL;  
00487   //  point->label_level=NULL;
00488   
00489 }
00490   */
00491 
00492   type = Give_type(I);
00493 
00494   if(type==multigrid || type==parallel_p) {
00495     if(point->label_level==NULL) {
00496       // initialize label_level since max_num_label==1      
00497       if(developer_version) {
00498         if(Give_max_num_label()!=1) 
00499           cout << " Error 2 in Set_label_bo_mg! ";
00500       }
00501       point->label_level = new int[1];
00502       point->label_level[0] = 0;
00503     }
00504 
00505     if(label==true) {
00506       // Dirichlet auf level
00507       if(point->Give_label(num) < level)
00508         point->Put_label_level(num,level);
00509     }
00510     // in case of label==false (Neumann) nothing has to be done
00511     //    }
00512 
00513     //  if(type==parallel_p) {
00514     if(label==true)
00515       point->send_label = Dirichlet_send;
00516     else 
00517       point->send_label = Neumann_send;
00518   }
00519 }
00520 
00521 bool Grid_base::Give_label_bo_D(Index3D I, int num, int level) const {
00522   Point_hashtable1* point;
00523   
00524   point=hashtable1_point(I);
00525   if(developer_version) {
00526     if(point==NULL) {
00527       cout << " Error in Give_label_bo(Index3D Ind, int num)" << endl;
00528     }
00529   }
00530   return (point->Give_label(num)>=level);
00531 }
00532 
00533 
00534 // 3. Parallelization
00536 
00537 
00538 
00539 
00540 bool Grid_base::Send_label_in_direction_A(int d, 
00541                                           Point_hashtable1* poi,
00542                                           int t, 
00543                                           Index3D my_lev_index,
00544                                           Index3D next_index) {
00545   Index3D I;
00546   Index3D nextS;
00547 
00548   /*
00549   I = poi->ind;
00550   if(poi->send_label!=dont_send && t>=I.Tiefe()) {
00551     if(I << next_index)
00552       return true;
00553   }
00554   */
00555  I = poi->ind;
00556   if(poi->send_label!=dont_send && t>=I.Tiefe()) {
00557     if(I << next_index)
00558       return true;
00559     
00560     if(d<6) {
00561       nextS = I.next((dir_3D)d,t);
00562     }
00563     else if(d<18) {
00564       nextS = I.next((Edges_cell)(d-6),t);
00565     }
00566     else
00567       nextS = I.next((dir_sons)(d-18),t);
00568     
00569     if(nextS << next_index)
00570       return true;
00571   }
00572 
00573   return false;
00574 }
00575 
00576 bool Grid_base::Send_label_in_direction_B(int d, 
00577                                           Point_hashtable1* poi,
00578                                           int t, 
00579                                           Index3D my_lev_index,
00580                                           Index3D next_index) {
00581   Index3D I;
00582   Index3D nextS;
00583 
00584   I = poi->ind;
00585   if(poi->send_label!=dont_send && t>=I.Tiefe()) {
00586     if(I << next_index)
00587       return true;
00588     
00589     if(d<6) {
00590       nextS = I.next((dir_3D)d,t);
00591     }
00592     else if(d<18) {
00593       nextS = I.next((Edges_cell)(d-6),t);
00594     }
00595     else
00596       nextS = I.next((dir_sons)(d-18),t);
00597     
00598     if(nextS << next_index)
00599       return true;
00600   }
00601 
00602   return false;
00603 }
00604 
00605 
00606 
00607 // level ist hier Punkte level
00608 void Grid_base::Send_label_parallel(int level, int num_label) {
00609   int d, level_label;
00610 
00611   Index3D I, next_index;
00612   int number_bo_send, number_bo_rec, num_message, num;
00613   Index3D my_lev_index;
00614 
00615   int* send_info;
00616   int* receive_info;
00617 
00618   MPI_Request req[2];
00619   MPI_Status status[2];
00620 
00621   int rank_source, rank_destination;
00622 
00623   // Part A: from fine to coarse
00624   //-----------------------------
00625 
00626   /*
00627   cout << " Send_label_parallel! level: " 
00628        << level 
00629        << " maxlevel: " << max_level 
00630        << " my_rank: " << my_rank
00631        << endl;
00632   */
00633 
00634   // -1. set my_lev_index 
00635   my_lev_index = Give_my_level_index(level+2);
00636   //my_lev_index = Give_my_level_index(level+1);
00637 
00638   for(d=0;d<26;++d) { 
00639     
00640     // 0. part: basic infos about send and recieve
00641     if(d<6) {
00642       rank_source      = give_next_rank_source((dir_3D)d,level+1);
00643       rank_destination = give_next_rank_destination((dir_3D)d,level+1);  
00644     }
00645     else if(d<18) {
00646       rank_source      = give_next_rank_source((Edges_cell)(d-6),level+1);
00647       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level+1);  
00648     }
00649     else {
00650       rank_source      = give_next_rank_source((dir_sons)(d-18),level+1);
00651       rank_destination = give_next_rank_destination((dir_sons)(d-18),level+1);
00652     }
00653     
00654 
00655     /*    
00656     if(d<6) {
00657       rank_source      = give_next_rank_source((dir_3D)d,level);
00658       rank_destination = give_next_rank_destination((dir_3D)d,level);  
00659     }
00660     else if(d<18) {
00661       rank_source      = give_next_rank_source((Edges_cell)(d-6),level);
00662       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level);  
00663     }
00664     else {
00665       rank_source      = give_next_rank_source((dir_sons)(d-18),level);
00666       rank_destination = give_next_rank_destination((dir_sons)(d-18),level);
00667     }
00668     */
00669 
00670     if(developer_version) {
00671       if(rank_source>p || rank_source<-1) {
00672         cout << " Problem in Send_label_parallel! " 
00673              << " level: " << level 
00674              << " d: " << d << endl;
00675       }
00676     }
00677     
00678 
00679     // 1. index of receiving next processor and correct rank's
00680     //next_index = Give_next_on_level(d,my_lev_index);
00681     next_index = Give_level_index_of(rank_destination,level+1);
00682 
00683     // periodisch neu
00684     if(next_index == Index3D(3,3,3))
00685       rank_destination = -1;
00686  
00687     if(Give_my_coarsest_level()>=level+1)
00688       rank_source      = -1;
00689 
00690     // 2. part: count number of points to be send
00691     number_bo_send = 0;
00692     if(rank_destination != -1) {
00693       iterate_hash1 {
00694 
00695         /*
00696         //      if(point1->ind==Index3D(0,2,0)) {
00697         if(point1->ind==Index3D(37,17,33)) {
00698           //      if(my_rank==17)
00699             cout << " destination A: " << rank_destination 
00700                  << " level: " << level 
00701                  << " wie: " <<
00702               Send_label_in_direction_A(d, point1,level,
00703                                       my_lev_index,next_index)
00704                  << " label: " <<  (point1->send_label!=dont_send)
00705                  << " incl: " << (point1->ind << next_index)
00706                  << " x: " << next_index.I_x().get()
00707                  << " y: " << next_index.I_y().get()
00708                  << " z: " << next_index.I_z().get()
00709                  << " my_rank: " << my_rank
00710                  << endl;
00711         }
00712         */
00713         
00714 
00715         if(Send_label_in_direction_A(d, point1,level,
00716                                      my_lev_index,next_index)) {
00717           ++number_bo_send;
00718         }
00719       }
00720     }
00721     
00722     // 3. part: send informations about number_bo_send
00723     number_bo_rec = 0;
00724     num_message = 0;
00725     if(rank_source != -1) {
00726       MPI_Irecv(&number_bo_rec,1,MPI_INT,rank_source,100+d,comm,
00727                 &req[num_message]);
00728       ++num_message;
00729     }
00730     if(rank_destination != -1) {
00731       MPI_Isend(&number_bo_send,1,MPI_INT,rank_destination,100+d,comm,
00732                 &req[num_message]);
00733       ++num_message;
00734     }
00735     MPI_Waitall(num_message,req,status);
00736     
00737     // 4. part: build informations about points to be send
00738     // How informations are stored:
00739     // Index of cell       
00740     // ind_x, ind_y, ind_z, level(Dirichlet)     >   4 int
00741     send_info = Give_buffer_send_info(4*number_bo_send);
00742     receive_info = Give_buffer_receive_info(4*number_bo_rec);
00743     
00744     num = 0; 
00745     if(rank_destination != -1) {
00746       iterate_hash1 {
00747         if(Send_label_in_direction_A(d, point1,level,
00748                                      my_lev_index,next_index)) {
00749           I = point1->ind;
00750           
00751           send_info[num*4]   = I.I_x().get();
00752           send_info[num*4+1] = I.I_y().get();
00753           send_info[num*4+2] = I.I_z().get();
00754           
00755           if(point1->send_label==Dirichlet_send)
00756             send_info[num*4+3] = level;
00757           else 
00758             send_info[num*4+3] = 0;
00759           num++;
00760         }
00761       }
00762     }
00763     
00764     // 4. part: send points
00765     num_message = 0;
00766     if(rank_source != -1 && number_bo_rec>0) {
00767       MPI_Irecv(receive_info,4*number_bo_rec,MPI_INT,rank_source,120+d,comm,
00768                 &req[num_message]);
00769       ++num_message;
00770     }
00771     if(rank_destination != -1 && number_bo_send>0) {
00772       MPI_Isend(send_info,4*number_bo_send,MPI_INT,rank_destination,120+d,comm,
00773                 &req[num_message]);
00774       ++num_message;
00775     }
00776     // MPI_Barrier(comm);
00777     MPI_Waitall(num_message,req,status);
00778     
00779     // 5. part: set points
00780     for(num=0;num<number_bo_rec;++num) {
00781       I = Index3D(receive_info[num*4],
00782                   receive_info[num*4+1],
00783                   receive_info[num*4+2]);
00784       level_label = receive_info[num*4+3];
00785       if(level_label!=0)
00786         Set_label_bo_mg(I, num_label, true, level_label);
00787       else 
00788         Set_label_bo_mg(I, num_label, false, level_label);
00789 
00790 
00791 
00792 
00793       /*
00794         if(I==Index3D(37,17,33)) {
00795           cout << " get A: " << rank_destination 
00796                << " level: " << level 
00797                << " incl: " << level_label
00798                << " my_rank: " << my_rank
00799                  << endl;
00800         }
00801       */
00802 
00803 
00804 
00805 
00806 
00807     }
00808   }
00809 
00810 
00811   /*
00812   // Part B: ditribute on coarse
00813   //-----------------------------
00814 
00815   // -1. set my_lev_index 
00816   my_lev_index = Give_my_level_index(level+1);
00817 
00818   for(d=0;d<26;++d) { 
00819     
00820     // 0. part: basic infos about send and recieve
00821     if(d<6) {
00822       rank_source      = give_next_rank_source((dir_3D)d,level);
00823       rank_destination = give_next_rank_destination((dir_3D)d,level);  
00824     }
00825     else if(d<18) {
00826       rank_source      = give_next_rank_source((Edges_cell)(d-6),level);
00827       rank_destination = give_next_rank_destination((Edges_cell)(d-6),level);  
00828     }
00829     else {
00830       rank_source      = give_next_rank_source((dir_sons)(d-18),level);
00831       rank_destination = give_next_rank_destination((dir_sons)(d-18),level);
00832     }
00833 
00834     // 1. index of receiving next processor and correct rank's
00835     next_index = Give_next_on_level(d,my_lev_index);
00836 
00837     // 2. part: count number of points to be send
00838     number_bo_send = 0;
00839     if(rank_destination != -1) {
00840       iterate_hash1 {
00841         if(Send_label_in_direction_B(d, point1,level,
00842                                      my_lev_index,next_index)) {
00843           ++number_bo_send;
00844         }
00845       }
00846     }
00847     
00848     // 3. part: send informations about number_bo_send
00849     number_bo_rec = 0;
00850     num_message = 0;
00851     if(rank_source != -1) {
00852       MPI_Irecv(&number_bo_rec,1,MPI_INT,rank_source,100+d,comm,
00853                 &req[num_message]);
00854       ++num_message;
00855     }
00856     if(rank_destination != -1) {
00857       MPI_Isend(&number_bo_send,1,MPI_INT,rank_destination,100+d,comm,
00858                 &req[num_message]);
00859       ++num_message;
00860     }
00861     MPI_Waitall(num_message,req,status);
00862     
00863     // 4. part: build informations about points to be send
00864     // How informations are stored:
00865     // Index of cell       
00866     // ind_x, ind_y, ind_z, level(Dirichlet)     >   4 int
00867     send_info = Give_buffer_send_info(4*number_bo_send);
00868     receive_info = Give_buffer_receive_info(4*number_bo_rec);
00869     
00870     num = 0; 
00871     if(rank_destination != -1) {
00872       iterate_hash1 {
00873         if(Send_label_in_direction_B(d, point1,level,
00874                                      my_lev_index,next_index)) {
00875           I = point1->ind;
00876           
00877           send_info[num*4]   = I.I_x().get();
00878           send_info[num*4+1] = I.I_y().get();
00879           send_info[num*4+2] = I.I_z().get();
00880           
00881           if(point1->send_label==Dirichlet_send)
00882             send_info[num*4+3] = level;
00883           else 
00884             send_info[num*4+3] = 0;
00885           num++;
00886         }
00887       }
00888     }
00889     
00890     // 4. part: send points
00891     num_message = 0;
00892     if(rank_source != -1 && number_bo_rec>0) {
00893       MPI_Irecv(receive_info,4*number_bo_rec,MPI_INT,rank_source,120+d,comm,
00894                 &req[num_message]);
00895       ++num_message;
00896     }
00897     if(rank_destination != -1 && number_bo_send>0) {
00898       MPI_Isend(send_info,4*number_bo_send,MPI_INT,rank_destination,120+d,comm,
00899                 &req[num_message]);
00900       ++num_message;
00901     }
00902     // MPI_Barrier(comm);
00903     MPI_Waitall(num_message,req,status);
00904     
00905     // 5. part: set points
00906     for(num=0;num<number_bo_rec;++num) {
00907       I = Index3D(receive_info[num*4],
00908                   receive_info[num*4+1],
00909                   receive_info[num*4+2]);
00910       level_label = receive_info[num*4+3];
00911       if(level_label!=0)
00912         Set_label_bo_mg(I, num_label, true, level_label);
00913       else 
00914         Set_label_bo_mg(I, num_label, false, level_label);
00915     }
00916   }
00917   */
00918 }
00919 
00920 
00921 
00922 
00923 

Generated on Mon Jan 16 13:23:41 2006 for IPPL by  doxygen 1.4.6