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

src/expde/extemp/variabl2.h

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 // variable2.h
00020 //
00021 // ------------------------------------------------------------
00022 
00023 #ifndef VARIABLEZ_H_
00024 #define VARIABLEZ_H_
00025                 
00026 // ------------------------------------------------------------
00027 //  implementation of memberfunctions with are not inline
00028 // ------------------------------------------------------------
00029 
00030 
00031 
00032 
00033 template<class R>
00034 void Variable::operator=(const DExprResVarP<R> &vr) {
00035   bool run;
00036   int color;
00037   int run_bo;
00038   int lev;
00039   double *vars;
00040   
00041   Test_init();  // Grid storage initialization
00042 
00043   if(grid->I_am_active()) {
00044     lev = vr.Level();
00045     
00046     if(lev < 0)          { cout << "\n Level too coarse! "   << endl; }
00047     else if (lev > Max_Level())  { cout << "\n Level too fine! " << endl; }
00048     run = vr.run_interior();
00049     
00050     // for update type
00051     if(parallel_version) {
00052       grid->Set_type_of_update(number_variable,lev,no_update);
00053     }
00054     
00055     if(run && number_variable != vr.Number_variable()) {
00056       // innere Punkte
00057       for(color=0;color<8;++color)
00058         for(iter_i=grid->Start_P_interior(lev,color);
00059             iter_i!=NULL;iter_i=iter_i->Next()) {
00060           vars = iter_i->varM(grid,lev);
00061           vars[number_variable] = vars[vr.Number_variable()];
00062         }
00063     }
00064     
00065     Active_interior(run);
00066     run = vr.run_nearb();
00067     if(run && number_variable != vr.Number_variable()) {
00068       // Randzellfreiheitsgrade
00069       for(iter_cf=grid->Start_P_cellpoi(lev);
00070           iter_cf!=NULL;iter_cf=iter_cf->Next()) {
00071         vars = iter_cf->varM(grid);
00072         vars[number_variable] = vars[vr.Number_variable()];
00073       }
00074       // randnahe Punkte
00075       for(color=0;color<8;++color)
00076         for(iter_n=grid->Start_P_nearb(lev,color);
00077             iter_n!=NULL;iter_n=iter_n->Next()) {
00078           vars = iter_n->varM(grid,lev);
00079           vars[number_variable] = vars[vr.Number_variable()];
00080         }
00081     }
00082     Active_nearb(run);
00083     run_bo = vr.run_boundary();
00084     
00085     if(number_variable != vr.Number_variable()) {
00086       if(run_bo==true) {
00087         // Randpunkte
00088         for(color=0;color<8;++color) {
00089           for(iter_b=grid->Start_P_Bo2p(lev,color);
00090               iter_b!=NULL;iter_b=iter_b->Next()) {
00091             vars = iter_b->varM(grid);
00092             vars[number_variable] = vars[vr.Number_variable()];
00093           }
00094         }
00095         // exterior grobe Randpunkte
00096         for(color=0;color<8;++color) {
00097           for(iter_n=grid->Start_P_exteri(lev,color);
00098               iter_n!=NULL;iter_n=iter_n->Next()) {
00099             vars = iter_n->varM(grid,lev);
00100             vars[number_variable] = vars[vr.Number_variable()];
00101           }
00102         }
00103       }
00104       else if(run_bo<false) {
00105         // Randpunkte
00106         for(color=0;color<8;++color)
00107           for(iter_b=grid->Start_P_Bo2p(lev,color);
00108               iter_b!=NULL;iter_b=iter_b->Next()) {
00109             if(iter_b->Give_Label(-run_bo,grid)) {      
00110               vars = iter_b->varM(grid);
00111               vars[number_variable] = vars[vr.Number_variable()];
00112             }
00113           }
00114       }
00115     }
00116     
00117     Active_boundary(run_bo);
00118     Active_Level(lev);   
00119   }
00120 }
00121 
00122 
00123 
00124 // for assign; usual
00125 template<class A>
00126 void Variable::operator=(const DExpr<A>& a_) {
00127   int lev;
00128   int run_bo;
00129   double h_mesh;
00130 
00131 
00132   // 1. Step: Same preparations
00133   if(a_.Dominant_poi()==anti_dominant) run_bo = r_boundary;
00134   else                                 run_bo = a_.run_boundary();
00135 
00136   if(grid->I_am_active()) {
00137     a_.Put_grid_rbo(grid,run_bo);    // Give everybody the pointer to the grid
00138     Test_init();          // Grid storage initialization
00139   }
00140   // End preparations
00141 
00142   // 2. Step: level and meshsize
00143   if(a_.Dominant_lev()==anti_dominant) lev = level;
00144   else                                 lev = a_.Level();
00145   if(lev < 0)          { cout << "\n Level too coarse! "   << endl; }
00146   else if (lev > Max_Level())  { cout << "\n Level too fine! " << endl; }
00147   
00148   h_mesh = grid->H_mesh() / Zweipotenz(lev);
00149 
00150   Active_Level(lev);
00151   a_.Active_Sim_Level(lev);                    // For operator ==
00152 
00153 
00154   // 3. Step: activity of points
00155  if(grid->I_am_active()) {
00156    // do this only if active
00157   Run_type_object run_obj;
00158   
00159   if(a_.Dominant_poi()==anti_dominant) {
00160     run_obj.r_interior = r_interior;
00161     run_obj.r_nearb    = r_nearb;
00162     run_obj.r_boundary = r_boundary;
00163   }
00164   else {
00165     run_obj.r_interior = a_.run_interior();
00166     run_obj.r_nearb    = a_.run_nearb();
00167     run_obj.r_boundary = a_.run_boundary();
00168   }
00169   Active_interior(run_obj.r_interior);
00170   a_.Active_Sim_interior(run_obj.r_interior);  // For operator ==
00171   Active_nearb(run_obj.r_nearb);
00172   a_.Active_Sim_nearb(run_obj.r_nearb);        // For operator ==
00173   Active_boundary(run_obj.r_boundary);
00174   a_.Active_Sim_boundary(run_obj.r_boundary);  // For operator ==
00175  
00176 
00177   // 4.Step: Evaluation_Parallelization
00178   Evaluation_Parallelization_object* ev_par_obj;
00179   evaluation_typ et;
00180   
00181   // 4.1 prepare
00182   ev_par_obj = grid->Give_eval_par();
00183   // -> parallel
00184   if(parallel_version) {
00185     ev_par_obj->StartA_evaluation(run_obj,lev,a_.Sice_stencil(),
00186                                   a_.ops_interior());
00187     a_.Add_variables_for_parallel(ev_par_obj);  // this adds the variables
00188                                                 // and tells if communication
00189                                                 // for restriction and
00190                                                 // prolongation is necessary
00191     bool GS_expr;
00192     GS_expr = a_.GS_type(number_variable,0);
00193     ev_par_obj->StartB_evaluation(number_variable,GS_expr);
00194     a_.Active_Sim_update(ev_par_obj,
00195                          ev_par_obj->Give_left_level(), // For operator ==
00196                          ev_par_obj->Give_left_update()); // For operator ==
00197 
00198   }
00199   // -> seriel
00200   else {
00201     ev_par_obj->Start_evaluation_seriel(run_obj,lev,a_.Sice_stencil(),
00202                                         a_.ops_interior());
00203   }
00204   // 4.2 iterate
00205   for(et=ev_par_obj->Give_next_evaluation_typ();
00206       et!=stop_evaluation;
00207       et=ev_par_obj->Give_next_evaluation_typ()) {
00208     switch(et) {
00209     case interior_15_evaluation_fine: {
00210       evaluate_interior_15_evaluation_fine(a_, ev_par_obj, h_mesh, lev);
00211       break;
00212     }
00213     case interior_17_evaluation_fine: {
00214       evaluate_interior_17_evaluation_fine(a_, ev_par_obj, h_mesh, lev);
00215       break;
00216     }
00217     case interior_25_evaluation_fine: {
00218       evaluate_interior_25_evaluation_fine(a_, ev_par_obj, h_mesh, lev);
00219       break;
00220     }
00221     case interior_1_evaluation_fine: {
00222       evaluate_interior_1_evaluation_fine(a_, ev_par_obj, h_mesh, lev);
00223       break;
00224     }
00225     case interior_15_evaluation_coarse: {
00226       for(iter_i=ev_par_obj->Start_interior_pointer();
00227           iter_i!=NULL;iter_i=iter_i->Next()) {
00228          iter_i->varM(grid,lev)[number_variable] = 
00229             a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00230                                     iter_i->varN(grid,lev),
00231                                     iter_i->varW(grid,lev),
00232                                     iter_i->varM(grid,lev),
00233                                     iter_i->varE(grid,lev),
00234                                     iter_i->varS(grid,lev),
00235                                     iter_i->varT(grid,lev),
00236                                     iter_i->varD(grid,lev),
00237                                     iter_i->varND(grid,lev),
00238                                     iter_i->varWN(grid,lev),
00239                                     iter_i->varWT(grid,lev),
00240                                     iter_i->varED(grid,lev),
00241                                     iter_i->varST(grid,lev),
00242                                     iter_i->varES(grid,lev),
00243                                     NULL, NULL,
00244                                     iter_i->varEST(grid,lev),
00245                                     iter_i->varWND(grid,lev),
00246                                     NULL, NULL, NULL, NULL,
00247                                     NULL, NULL, NULL, NULL);
00248       }
00249       break;
00250     }
00251     case interior_25_evaluation_coarse: {
00252       for(iter_i=ev_par_obj->Start_interior_pointer();
00253           iter_i!=NULL;iter_i=iter_i->Next()) {
00254         iter_i->varM(grid,lev)[number_variable] = 
00255           a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00256                                   iter_i->varN(grid,lev),
00257                                   iter_i->varW(grid,lev),
00258                                   iter_i->varM(grid,lev),
00259                                   iter_i->varE(grid,lev),
00260                                   iter_i->varS(grid,lev),
00261                                   iter_i->varT(grid,lev),
00262                                   iter_i->varD(grid,lev),
00263                                   iter_i->varND(grid,lev),
00264                                   iter_i->varWN(grid,lev),
00265                                   iter_i->varWT(grid,lev),
00266                                   iter_i->varED(grid,lev),
00267                                   iter_i->varST(grid,lev),
00268                                   iter_i->varES(grid,lev),
00269                                   iter_i->varWD(grid,lev),
00270                                   iter_i->varET(grid,lev),
00271                                   iter_i->varEST(grid,lev),
00272                                   iter_i->varWND(grid,lev),
00273                                   iter_i->cell_varWSD(grid,lev),
00274                                   iter_i->cell_varWST(grid,lev),
00275                                   iter_i->cell_varWND(grid,lev),
00276                                   iter_i->cell_varWNT(grid,lev),
00277                                   iter_i->cell_varESD(grid,lev),
00278                                   iter_i->cell_varEST(grid,lev),
00279                                   iter_i->cell_varEND(grid,lev),
00280                                   iter_i->cell_varENT(grid,lev));
00281 
00282       }
00283       break;
00284     }
00285     case interior_17_evaluation_coarse: {
00286       for(iter_i=ev_par_obj->Start_interior_pointer();
00287           iter_i!=NULL;iter_i=iter_i->Next()) {
00288         iter_i->varM(grid,lev)[number_variable] = 
00289           a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00290                                   iter_i->varN(grid,lev),
00291                                   iter_i->varW(grid,lev),
00292                                   iter_i->varM(grid,lev),
00293                                   iter_i->varE(grid,lev),
00294                                   iter_i->varS(grid,lev),
00295                                   iter_i->varT(grid,lev),
00296                                   iter_i->varD(grid,lev),
00297                                   iter_i->varND(grid,lev),
00298                                   iter_i->varWN(grid,lev),
00299                                   iter_i->varWT(grid,lev),
00300                                   iter_i->varED(grid,lev),
00301                                   iter_i->varST(grid,lev),
00302                                   iter_i->varES(grid,lev),
00303                                   iter_i->varWD(grid,lev),
00304                                   iter_i->varET(grid,lev),
00305                                   iter_i->varEST(grid,lev),
00306                                   iter_i->varWND(grid,lev),
00307                                   NULL, NULL, NULL, NULL,
00308                                   NULL, NULL, NULL, NULL);
00309       }
00310       break;
00311     }
00312     case interior_1_evaluation_coarse: {
00313       for(iter_i=ev_par_obj->Start_interior_pointer();
00314           iter_i!=NULL;iter_i=iter_i->Next()) {
00315         iter_i->varM(grid,lev)[number_variable] = 
00316           a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00317                                   NULL, NULL,
00318                                   iter_i->varM(grid,lev),
00319                                   NULL, NULL, NULL, NULL,
00320                                   NULL, NULL, NULL, NULL,
00321                                   NULL, NULL,
00322                                   NULL, NULL, NULL, NULL,
00323                                   NULL, NULL, NULL, NULL,
00324                                   NULL, NULL, NULL, NULL);
00325       }
00326       break;
00327     }
00328     case cellf_large_evaluation: {
00329       for(iter_cf=ev_par_obj->Start_cellpoi_pointer();
00330           iter_cf!=NULL;iter_cf=iter_cf->Next()) {
00331         iter_cf->varM(grid)[number_variable] = 
00332           a_.Give_cellpoi(iter_cf,grid,grid->Give_Bo_cell(iter_cf->Ind()));
00333       }
00334       break;
00335     }
00336     case cellf_1_evaluation: {
00337       for(iter_cf=ev_par_obj->Start_cellpoi_pointer();
00338           iter_cf!=NULL;iter_cf=iter_cf->Next()) {
00339         iter_cf->varM(grid)[number_variable] = 
00340           a_.Give_cellpoi(iter_cf,grid,NULL);
00341       }
00342       break;
00343     }
00344     case nearb_15_evaluation: {
00345       nearb_ablage = grid->Give_Nearb_Ablage();
00346       for(iter_n=ev_par_obj->Start_nearb_pointer();
00347           iter_n!=NULL;iter_n=iter_n->Next()) {
00348         nearb_ablage->Initialize(grid,iter_n,lev);
00349         iter_n->varM(grid,lev)[number_variable] = 
00350           a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00351       }
00352       break;
00353     }
00354     case nearb_17_evaluation: {
00355       nearb_ablage = grid->Give_Nearb_Ablage();
00356       for(iter_n=ev_par_obj->Start_nearb_pointer();
00357           iter_n!=NULL;iter_n=iter_n->Next()) {
00358         nearb_ablage->Initialize_17(grid,iter_n,lev);
00359         iter_n->varM(grid,lev)[number_variable] = 
00360           a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00361       }
00362       break;
00363     }
00364     case nearb_1_evaluation: {
00365       for(iter_n=ev_par_obj->Start_nearb_pointer();
00366           iter_n!=NULL;iter_n=iter_n->Next()) {
00367         iter_n->varM(grid,lev)[number_variable] = 
00368           a_.Give_nearb(iter_n,grid,h_mesh,lev,NULL);
00369       }
00370       break;
00371     }
00372     case boundary_evaluation_all: {
00373       for(iter_b=ev_par_obj->Start_Bo2p_pointer();
00374           iter_b!=NULL;iter_b=iter_b->Next()) {
00375         iter_b->varM(grid)[number_variable] = 
00376           a_.Give_Bo2p(iter_b,grid,lev);
00377       }
00378       break;
00379     }
00380     case boundary_evaluation_subset: {
00381       for(iter_b=ev_par_obj->Start_Bo2p_pointer();
00382           iter_b!=NULL;iter_b=iter_b->Next()) {
00383         if(iter_b->Give_Label(-run_bo,grid))
00384           iter_b->varM(grid)[number_variable] = 
00385             a_.Give_Bo2p(iter_b,grid,lev);
00386       }
00387       break;
00388     }
00389     case exteri_1_evaluation_subset: {
00390       for(iter_n=ev_par_obj->Start_nearb_pointer();
00391           iter_n!=NULL;iter_n=iter_n->Next()) {
00392         if(iter_n->Give_Label(-run_bo,lev,grid))
00393           iter_n->varM(grid,lev)[number_variable] = 
00394             a_.Give_nearb(iter_n,grid,h_mesh,lev,NULL);
00395       }
00396       break;
00397     }
00398     case exteri_17_evaluation_subset: {
00399       nearb_ablage = grid->Give_Nearb_Ablage();
00400       for(iter_n=ev_par_obj->Start_nearb_pointer();
00401           iter_n!=NULL;iter_n=iter_n->Next()) {
00402         if(iter_n->Give_Label(-run_bo,lev,grid)) {
00403           nearb_ablage->Initialize_17(grid,iter_n,lev);
00404           iter_n->varM(grid,lev)[number_variable] = 
00405             a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00406         }
00407       }
00408       break;
00409     }
00410     case exteri_15_evaluation_subset: {
00411       nearb_ablage = grid->Give_Nearb_Ablage();
00412       for(iter_n=ev_par_obj->Start_nearb_pointer();
00413           iter_n!=NULL;iter_n=iter_n->Next()) {
00414         if(iter_n->Give_Label(-run_bo,lev,grid)) {
00415           nearb_ablage->Initialize(grid,iter_n,lev);
00416           iter_n->varM(grid,lev)[number_variable] = 
00417             a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00418         }
00419       }
00420       break;
00421     }
00422     default: {
00423       break;
00424     }
00425     }
00426   }
00427   
00428   // 5.Step: clean storage
00429   grid->Give_Nearb_Ablage()->Put_leer();
00430 
00431   a_.clean();
00432  }
00433 }
00434 
00435 
00436 
00437 
00438 // for assign; usual
00439 template<class A>
00440 void Variable::evaluate_interior_15_evaluation_fine(const DExpr<A>& a_,
00441                         Evaluation_Parallelization_object* ev_par_obj,
00442                                                     double h_mesh, int lev)  {
00443   for(iter_i=ev_par_obj->Start_interior_pointer();
00444       iter_i!=NULL;iter_i=iter_i->Next()) {
00445     iter_i->varM(grid,lev)[number_variable] = 
00446       a_.Give_interior(iter_i,grid,h_mesh,lev,
00447                        iter_i->varN(grid,lev),
00448                        iter_i->varW(grid,lev),
00449                        iter_i->varM(grid,lev),
00450                        iter_i->varE(grid,lev),
00451                        iter_i->varS(grid,lev),
00452                        iter_i->varT(grid,lev),
00453                        iter_i->varD(grid,lev),
00454                        iter_i->varND(grid,lev),
00455                        iter_i->varWN(grid,lev),
00456                        iter_i->varWT(grid,lev),
00457                        iter_i->varED(grid,lev),
00458                        iter_i->varST(grid,lev),
00459                        iter_i->varES(grid,lev),
00460                        NULL, NULL,
00461                        iter_i->varEST(grid,lev),
00462                        iter_i->varWND(grid,lev),
00463                        NULL, NULL, NULL, NULL,
00464                        NULL, NULL, NULL, NULL);
00465   }
00466 }
00467 
00468 
00469 template<class A>
00470 void Variable::evaluate_interior_17_evaluation_fine(const DExpr<A>& a_,
00471                         Evaluation_Parallelization_object* ev_par_obj,
00472                                                     double h_mesh, int lev)  {
00473   for(iter_i=ev_par_obj->Start_interior_pointer();
00474       iter_i!=NULL;iter_i=iter_i->Next()) {
00475     iter_i->varM(grid,lev)[number_variable] = 
00476       a_.Give_interior(iter_i,grid,h_mesh,lev,
00477                        iter_i->varN(grid,lev),
00478                        iter_i->varW(grid,lev),
00479                        iter_i->varM(grid,lev),
00480                        iter_i->varE(grid,lev),
00481                        iter_i->varS(grid,lev),
00482                        iter_i->varT(grid,lev),
00483                        iter_i->varD(grid,lev),
00484                        iter_i->varND(grid,lev),
00485                        iter_i->varWN(grid,lev),
00486                        iter_i->varWT(grid,lev),
00487                        iter_i->varED(grid,lev),
00488                        iter_i->varST(grid,lev),
00489                        iter_i->varES(grid,lev),
00490                        iter_i->varWD(grid,lev),
00491                        iter_i->varET(grid,lev),
00492                        iter_i->varEST(grid,lev),
00493                        iter_i->varWND(grid,lev),
00494                        NULL, NULL, NULL, NULL,
00495                        NULL, NULL, NULL, NULL);
00496   }
00497 }
00498 
00499 
00500 
00501 
00502 template<class A>
00503 void Variable::evaluate_interior_25_evaluation_fine(const DExpr<A>& a_,
00504                         Evaluation_Parallelization_object* ev_par_obj,
00505                                                     double h_mesh, int lev)  {
00506   for(iter_i=ev_par_obj->Start_interior_pointer();
00507       iter_i!=NULL;iter_i=iter_i->Next()) {
00508     iter_i->varM(grid,lev)[number_variable] = 
00509       a_.Give_interior(iter_i,grid,h_mesh,lev,
00510                        iter_i->varN(grid,lev),
00511                        iter_i->varW(grid,lev),
00512                        iter_i->varM(grid,lev),
00513                        iter_i->varE(grid,lev),
00514                        iter_i->varS(grid,lev),
00515                        iter_i->varT(grid,lev),
00516                        iter_i->varD(grid,lev),
00517                        iter_i->varND(grid,lev),
00518                        iter_i->varWN(grid,lev),
00519                        iter_i->varWT(grid,lev),
00520                        iter_i->varED(grid,lev),
00521                        iter_i->varST(grid,lev),
00522                        iter_i->varES(grid,lev),
00523                        iter_i->varWD(grid,lev),
00524                        iter_i->varET(grid,lev),
00525                        iter_i->varEST(grid,lev),
00526                        iter_i->varWND(grid,lev),
00527                        iter_i->cell_varWSD(grid,lev),
00528                        iter_i->cell_varWST(grid,lev),
00529                        iter_i->cell_varWND(grid,lev),
00530                        iter_i->cell_varWNT(grid,lev),
00531                        iter_i->cell_varESD(grid,lev),
00532                        iter_i->cell_varEST(grid,lev),
00533                        iter_i->cell_varEND(grid,lev),
00534                        iter_i->cell_varENT(grid,lev));
00535   }
00536 }
00537 
00538 
00539 
00540 
00541 
00542 template<class A>
00543 void Variable::evaluate_interior_1_evaluation_fine(const DExpr<A>& a_,
00544                         Evaluation_Parallelization_object* ev_par_obj,
00545                                                    double h_mesh, int lev)  {
00546   for(iter_i=ev_par_obj->Start_interior_pointer();
00547       iter_i!=NULL;iter_i=iter_i->Next()) {
00548     iter_i->varM(grid,lev)[number_variable] = 
00549       a_.Give_interior(iter_i,grid,h_mesh,lev,
00550                        NULL, NULL,
00551                        iter_i->varM(grid,lev),
00552                        NULL, NULL, NULL, NULL,
00553                        NULL, NULL, NULL, NULL,
00554                        NULL, NULL,
00555                        NULL, NULL, NULL, NULL,
00556                        NULL, NULL, NULL, NULL,
00557                        NULL, NULL, NULL, NULL);
00558   }
00559 }
00560 
00561 
00562 
00563 
00564 
00566 
00567 
00568 
00569 // for assign; restict domain
00570 template<class A>
00571 void Variable::operator=(const DExprResDomain<A>& a_) {
00572   int lev;
00573   double h_mesh;
00574 
00575   // Same preparations
00576   a_.Put_grid_rbo(grid,0);    // Give everybody the pointer to the grid
00577 
00578   Test_init();          // Grid storage initialization
00579   // End preparations
00580 
00581   if(a_.Dominant_lev()==anti_dominant) lev = level;
00582   else                                 lev = a_.Level();
00583   if(lev < 0)          { cout << "\n Level too coarse! "   << endl; }
00584   else if (lev > Max_Level())  { cout << "\n Level too fine! " << endl; }
00585   
00586   h_mesh = grid->H_mesh() / Zweipotenz(lev);
00587 
00588   evaluate_interior_dom(a_,lev,h_mesh);
00589   evaluate_nearb_dom(a_,lev,h_mesh);
00590   evaluate_boundary_dom(a_,lev);
00591 
00592   // clean storage
00593   a_.clean();
00594 }
00595 
00596 
00597 template<class A>
00598 void Variable::evaluate_interior_dom(const DExprResDomain<A>& a_, 
00599                                      int lev, double h_mesh) {
00600   // interior
00601   //----------
00602   grid->Start_calc_time_interior();  
00603   // ON  THE FINEST GRID
00604   if (lev == Max_Level()) {
00605     // number of operations
00606     grid->Add_operations(a_.ops_interior(),lev);
00607     // innere Punkte
00608     if(a_.Sice_stencil()==15) {  // 15-Stern
00609       for(iter_i=grid->Start_P_interior(lev);
00610           iter_i!=NULL;iter_i=iter_i->Next()) {
00611         if(a_.point_in_domain(iter_i->coordinateX(grid),
00612                               iter_i->coordinateY(grid),
00613                               iter_i->coordinateZ(grid)))
00614           iter_i->varM(grid,lev)[number_variable] = 
00615             a_.Give_interior(iter_i,grid,h_mesh,lev,
00616                              iter_i->varN(grid,lev),
00617                              iter_i->varW(grid,lev),
00618                              iter_i->varM(grid,lev),
00619                              iter_i->varE(grid,lev),
00620                              iter_i->varS(grid,lev),
00621                              iter_i->varT(grid,lev),
00622                              iter_i->varD(grid,lev),
00623                              iter_i->varND(grid,lev),
00624                              iter_i->varWN(grid,lev),
00625                              iter_i->varWT(grid,lev),
00626                              iter_i->varED(grid,lev),
00627                              iter_i->varST(grid,lev),
00628                              iter_i->varES(grid,lev),
00629                              NULL, NULL,
00630                              iter_i->varEST(grid,lev),
00631                              iter_i->varWND(grid,lev));
00632       }
00633     }
00634     else if(a_.Sice_stencil()==17) {  // 17-Stern
00635       for(iter_i=grid->Start_P_interior(lev);
00636           iter_i!=NULL;iter_i=iter_i->Next()) {
00637         if(a_.point_in_domain(iter_i->coordinateX(grid),
00638                               iter_i->coordinateY(grid),
00639                               iter_i->coordinateZ(grid)))
00640           iter_i->varM(grid,lev)[number_variable] = 
00641             a_.Give_interior(iter_i,grid,h_mesh,lev,
00642                              iter_i->varN(grid,lev),
00643                              iter_i->varW(grid,lev),
00644                              iter_i->varM(grid,lev),
00645                              iter_i->varE(grid,lev),
00646                              iter_i->varS(grid,lev),
00647                              iter_i->varT(grid,lev),
00648                              iter_i->varD(grid,lev),
00649                              iter_i->varND(grid,lev),
00650                              iter_i->varWN(grid,lev),
00651                              iter_i->varWT(grid,lev),
00652                              iter_i->varED(grid,lev),
00653                              iter_i->varST(grid,lev),
00654                              iter_i->varES(grid,lev),
00655                              iter_i->varWD(grid,lev),
00656                              iter_i->varET(grid,lev),
00657                              iter_i->varEST(grid,lev),
00658                              iter_i->varWND(grid,lev));
00659       }
00660     }
00661     else { // 1-Stern
00662       for(iter_i=grid->Start_P_interior(lev);
00663           iter_i!=NULL;iter_i=iter_i->Next()) {
00664         if(a_.point_in_domain(iter_i->coordinateX(grid),
00665                               iter_i->coordinateY(grid),
00666                               iter_i->coordinateZ(grid)))
00667           iter_i->varM(grid,lev)[number_variable] = 
00668             a_.Give_interior(iter_i,grid,h_mesh,lev,
00669                              NULL, NULL,
00670                              iter_i->varM(grid,lev),
00671                              NULL, NULL, NULL, NULL,
00672                              NULL, NULL, NULL, NULL,
00673                              NULL, NULL,
00674                              NULL, NULL, NULL, NULL);
00675       }
00676     }
00677   }
00678   else {
00679     // ON COARSER GRID
00680     // number of operations
00681     grid->Add_operations(a_.ops_interior(),lev);
00682     // innere Punkte
00683     h_mesh = grid->H_mesh() / Zweipotenz(lev);
00684     if(a_.Sice_stencil()==15) {  // 15-Stern
00685       for(iter_i=grid->Start_P_interior(lev);
00686           iter_i!=NULL;iter_i=iter_i->Next()) {
00687         if(a_.point_in_domain(iter_i->coordinateX(grid),
00688                               iter_i->coordinateY(grid),
00689                               iter_i->coordinateZ(grid)))
00690           iter_i->varM(grid,lev)[number_variable] = 
00691             a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00692                                     iter_i->varN(grid,lev),
00693                                     iter_i->varW(grid,lev),
00694                                     iter_i->varM(grid,lev),
00695                                     iter_i->varE(grid,lev),
00696                                     iter_i->varS(grid,lev),
00697                                     iter_i->varT(grid,lev),
00698                                     iter_i->varD(grid,lev),
00699                                     iter_i->varND(grid,lev),
00700                                     iter_i->varWN(grid,lev),
00701                                     iter_i->varWT(grid,lev),
00702                                     iter_i->varED(grid,lev),
00703                                     iter_i->varST(grid,lev),
00704                                     iter_i->varES(grid,lev),
00705                                     NULL, NULL,
00706                                     iter_i->varEST(grid,lev),
00707                                     iter_i->varWND(grid,lev));
00708       }
00709     }
00710     else if(a_.Sice_stencil()==17) {  // 17-Stern
00711       for(iter_i=grid->Start_P_interior(lev);
00712           iter_i!=NULL;iter_i=iter_i->Next()) {
00713         if(a_.point_in_domain(iter_i->coordinateX(grid),
00714                               iter_i->coordinateY(grid),
00715                               iter_i->coordinateZ(grid)))
00716           iter_i->varM(grid,lev)[number_variable] = 
00717             a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00718                                     iter_i->varN(grid,lev),
00719                                     iter_i->varW(grid,lev),
00720                                     iter_i->varM(grid,lev),
00721                                     iter_i->varE(grid,lev),
00722                                     iter_i->varS(grid,lev),
00723                                     iter_i->varT(grid,lev),
00724                                     iter_i->varD(grid,lev),
00725                                     iter_i->varND(grid,lev),
00726                                     iter_i->varWN(grid,lev),
00727                                     iter_i->varWT(grid,lev),
00728                                     iter_i->varED(grid,lev),
00729                                     iter_i->varST(grid,lev),
00730                                     iter_i->varES(grid,lev),
00731                                     iter_i->varWD(grid,lev),
00732                                     iter_i->varET(grid,lev),
00733                                     iter_i->varEST(grid,lev),
00734                                     iter_i->varWND(grid,lev));
00735       }
00736     }
00737     else { // 1-Stern
00738       for(iter_i=grid->Start_P_interior(lev);
00739           iter_i!=NULL;iter_i=iter_i->Next()) {
00740         if(a_.point_in_domain(iter_i->coordinateX(grid),
00741                               iter_i->coordinateY(grid),
00742                               iter_i->coordinateZ(grid)))
00743           iter_i->varM(grid,lev)[number_variable] = 
00744             a_.Give_interior_coarse(iter_i,grid,h_mesh,lev,
00745                                     NULL, NULL,
00746                                     iter_i->varM(grid,lev),
00747                                     NULL, NULL, NULL, NULL,
00748                                     NULL, NULL, NULL, NULL,
00749                                     NULL, NULL,
00750                                     NULL, NULL, NULL, NULL);
00751       }
00752     }
00753   }
00754   grid->Stop_calc_time_interior();
00755 }
00756 
00757 template<class A>
00758 void Variable::evaluate_nearb_dom(const DExprResDomain<A>& a_,
00759                                   int lev, double h_mesh) {
00760   // nearb
00761   //---------
00762   grid->Start_calc_time_nearbb();
00763   // Randzellfreiheitsgrade
00764   if(a_.Sice_stencil()>1) { // groesserer-Stern
00765     for(iter_cf=grid->Start_P_cellpoi(lev);
00766         iter_cf!=NULL;iter_cf=iter_cf->Next()) {
00767         if(a_.point_in_domain(iter_cf->coordinateX(grid),
00768                               iter_cf->coordinateY(grid),
00769                               iter_cf->coordinateZ(grid)))
00770           iter_cf->varM(grid)[number_variable] = 
00771             a_.Give_cellpoi(iter_cf,grid,grid->Give_Bo_cell(iter_cf->Ind()));
00772     }
00773   }
00774   else { // 1-Stern
00775     for(iter_cf=grid->Start_P_cellpoi(lev);
00776         iter_cf!=NULL;iter_cf=iter_cf->Next()) {
00777       if(a_.point_in_domain(iter_cf->coordinateX(grid),
00778                             iter_cf->coordinateY(grid),
00779                             iter_cf->coordinateZ(grid)))
00780         iter_cf->varM(grid)[number_variable] = 
00781           a_.Give_cellpoi(iter_cf,grid,NULL);
00782     }
00783   }
00784   
00785   // randnahe Punkte
00786   if(a_.Sice_stencil()>1) { // groesserer-Stern
00787     if(a_.Sice_stencil()==15) {
00788       nearb_ablage = grid->Give_Nearb_Ablage();
00789       for(iter_n=grid->Start_P_nearb(lev);
00790           iter_n!=NULL;iter_n=iter_n->Next()) {
00791         if(a_.point_in_domain(iter_n->coordinateX(grid),
00792                               iter_n->coordinateY(grid),
00793                               iter_n->coordinateZ(grid))) {
00794           nearb_ablage->Initialize(grid,iter_n,lev);
00795           iter_n->varM(grid,lev)[number_variable] = 
00796             a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00797         }
00798       }
00799     }
00800     else {
00801       nearb_ablage = grid->Give_Nearb_Ablage();
00802       for(iter_n=grid->Start_P_nearb(lev);
00803           iter_n!=NULL;iter_n=iter_n->Next()) {
00804         if(a_.point_in_domain(iter_n->coordinateX(grid),
00805                               iter_n->coordinateY(grid),
00806                               iter_n->coordinateZ(grid))) {
00807           nearb_ablage->Initialize_17(grid,iter_n,lev);
00808           iter_n->varM(grid,lev)[number_variable] = 
00809             a_.Give_nearb(iter_n,grid,h_mesh,lev,nearb_ablage);
00810         }
00811       }
00812     }
00813   }
00814   else { // 1-Stern
00815     for(iter_n=grid->Start_P_nearb(lev);
00816         iter_n!=NULL;iter_n=iter_n->Next()) {
00817       if(a_.point_in_domain(iter_n->coordinateX(grid),
00818                             iter_n->coordinateY(grid),
00819                             iter_n->coordinateZ(grid)))
00820         iter_n->varM(grid,lev)[number_variable] = 
00821           a_.Give_nearb(iter_n,grid,h_mesh,lev,NULL); 
00822     }
00823   }
00824   grid->Stop_calc_time_nearbb();
00825 }
00826 
00827 template<class A>
00828 void Variable::evaluate_boundary_dom(const DExprResDomain<A>& a_, 
00829                                      int lev) {
00830   // boundary
00831   //----------
00832   // Randpunkte
00833   grid->Start_calc_time_nearbb();
00834   if(a_.Sice_stencil()>1) { // groesserer-Stern
00835     for(iter_b=grid->Start_P_Bo2p(lev);
00836         iter_b!=NULL;iter_b=iter_b->Next()) {
00837       if(a_.point_in_domain(iter_b->coordinateX(grid),
00838                             iter_b->coordinateY(grid),
00839                             iter_b->coordinateZ(grid)))
00840         iter_b->varM(grid)[number_variable] = 
00841           a_.Give_Bo2p(iter_b,grid,lev);
00842     }
00843   }
00844   else { // 1-Stern
00845     for(iter_b=grid->Start_P_Bo2p(lev);
00846         iter_b!=NULL;iter_b=iter_b->Next()) {
00847       if(a_.point_in_domain(iter_b->coordinateX(grid),
00848                             iter_b->coordinateY(grid),
00849                             iter_b->coordinateZ(grid)))
00850         iter_b->varM(grid)[number_variable] = 
00851           a_.Give_Bo2p(iter_b,grid,lev);
00852     }   
00853   }
00854   grid->Stop_calc_time_nearbb();
00855 
00856   Active_Level(lev);   
00857   // For operator ==
00858   a_.Active_Sim_Level(lev);
00859 }
00860 
00861 
00862 
00863 #endif
00864 
00865             

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